Second (basic) tutorial¶
The H2 molecule, with convergence studies.¶
This tutorial aims at showing how to get converged values for the following physical properties:
- the bond length
- the atomisation energy
You will learn about the numerical quality of the calculations, then make convergence studies with respect to the number of planewaves and the size of the supercell, and finally consider the effect of the XC functional. The problems related to the use of different pseudopotential are not examined. You will also finish to read the abinit help file.
This tutorial should take about 1 hour.
Note
Supposing you made your own install of ABINIT, the input files to run the examples are in the ~abinit/tests/ directory where ~abinit is the absolute path of the abinit top-level directory. If you have NOT made your own install, ask your system administrator where to find the package, especially the executable and test files.
To execute the tutorials, create a working directory (Work*
) and
copy there the input files and the files file of the lesson. This will be explicitly mentioned in the first lessons,
that will tell you more about the files file (see also section 1.1).
The files file ending with _x (e.g. tbase1_x.files) must be edited every time you start to use a new input file.
Most of the tutorials do not rely on parallelism (except specific tutorials on parallelism). However you can run most of the tutorial examples in parallel, see the topic on parallelism.
In case you work on your own PC or workstation, to make things easier, we suggest you define some handy environment variables by executing the following lines in the terminal:
export ABI_HOME=Replace_with_the_absolute_path_to_the_abinit_top_level_dir export PATH=$ABI_HOME/src/98_main/:$PATH export ABI_TESTS=$ABI_HOME/tests/ export ABI_PSPDIR=$ABI_TESTS/Psps_for_tests/ # Pseudopotentials used in examples.
Examples in this tutorial use these shell variables: copy and paste the code snippets into the terminal (remember to set ABI_HOME first!). The ‘export PATH’ line adds the directory containing the executables to your PATH so that you can invoke the code by simply typing abinit in the terminal instead of providing the absolute path.
Summary of the previous tutorial¶
We studied the H_2 molecule in a big box. We used 10 Ha as cut-off energy, a 10x10x10 Bohr^3 supercell, the local-density approximation (as well as the local-spin-density approximation) in the Teter parametrization (ixc = 1, the default), and a pseudopotential from the Goedecker-Hutter-Teter table.
At this stage, we compared our results:
- bond length: 1.522 Bohr
- atomisation energy at that bond length: 0.1656 Ha = 4.506 eV
with the experimental data (as well as theoretical data using a much more accurate technique than DFT)
- bond length: 1.401 Bohr
- atomisation energy: 4.747 eV
The bond length is awful (nearly 10% off), and the atomisation energy is a bit too low, 5% off.
2 The convergence in ecut (I)¶
2.1.a Computing the bond length and corresponding atomisation energy in one run.
Before beginning, you might consider to work in a different subdirectory as for tutorial 1. Why not Work2?
Because we will compute many times the bond length and atomisation energy, it is worth to make a single input file that will do all the associated operations. You should try to use 2 datasets (try to combine $ABI_TESTS/tutorial/Input/tbase1_3.in with tbase1_5.in). Do not try to have the same position of the H atom as one of the H_2 atoms in the optimized geometry.
cd $ABI_TESTS/tutorial/Input mkdir Work2 cd Work2 cp ../tbase2_x.files . # You will need to edit this file. cp ../tbase2_1.in .
The input file tbase2_1.in is an example of file that will do the job,
# H2 molecule in a big box # # This file to optimize the H2 bond length, compute the associated total # energy, then to compute the total energy of the isolated H atom. # ndtset 2 #Definition of the unit cell and ecut, #for which one will have to make a convergence study acell 10 10 10 ecut 10 #First dataset : find the optimal bond length of H2, and associated total energy natom1 2 # There are two atoms ionmov1 2 # Use the modified Broyden algorithm ntime1 10 # Maximum number of Broyden "timesteps" tolmxf1 5.0d-4 # Stopping criterion for the geometry optimization : when # the residual forces are less than tolmxf, the Broyden # algorithm can stop xcart1 -0.7 0.0 0.0 # The starting values of the 0.7 0.0 0.0 # atomic coordinates toldff1 5.0d-5 # Will stop the SCF cycle when, twice in a row, # the difference between two consecutive evaluations of # forces differ by less than toldff (in Hartree/Bohr) nband1 1 # Just one band #Second dataset : get the total energy of the isolated atom natom2 1 # There is one atom nsppol2 2 # Spin-polarized calculation occopt2 2 # Allow occupation numbers to be set by hand nband2 1 1 # Number of bands for spin up and spin down occ2 1.0 0.0 # Occupation numbers for spin up state and spin down state. toldfe2 1.0d-6 # Will stop the SCF cycles when, twice in a row, # the difference between two consecutive evaluations # of total energy differ by less than toldfe (in Hartree) # This value is way too large for most realistic studies of materials xcart2 0.0 0.0 0.0 # The atom is located at the origin spinat2 0.0 0.0 1.0 # Initialisation of spin #rprim 1 0 0 0 1 0 0 0 1 # This line, defining orthogonal primitive vectors, # is commented, because it is precisely the default value of rprim #Definition of the atom types ntypat 1 # There is only one type of atom znucl 1 # The keyword "znucl" refers to the atomic number of the # possible type(s) of atom. The pseudopotential(s) # mentioned in the "files" file must correspond # to the type(s) of atom. Here, the only type is Hydrogen. #Definition of the atoms typat 1 1 # For the first dataset, both numbers will be read, # while for the second dataset, only one number will be read #Definition of the k-point grid kptopt 0 # Enter the k points manually nkpt 1 # Only one k point is needed for isolated system, # taken by default to be 0.0 0.0 0.0 #Definition of the SCF procedure nstep 10 # Maximal number of SCF cycles #toldfe is no more defined, as toldff is used above... diemac 2.0 # Although this is not mandatory, it is worth to # precondition the SCF cycle. The model dielectric # function used as the standard preconditioner # is described in the "dielng" input variable section. # Here, we follow the prescriptions for molecules # in a big box #%%<BEGIN TEST_INFO> #%% [setup] #%% executable = abinit #%% [files] #%% files_to_test = #%% tbase2_1.out, tolnlines= 3, tolabs= 3.000e-10, tolrel= 2.000e-08 #%% psp_files = 01h.pspgth #%% [paral_info] #%% max_nprocs = 1 #%% [extra_info] #%% authors = Unknown #%% keywords = #%% description = #%% H2 molecule in a big box #%% This file to optimize the H2 bond length, compute the associated total #%% energy, then to compute the total energy of the isolated H atom. #%%<END TEST_INFO>
while tbase2_1.out is an example of output file:
.Version 9.0.0 of ABINIT .(MPI version, prepared for a x86_64_linux_gnu9.2 computer) .Copyright (C) 1998-2020 ABINIT group . ABINIT comes with ABSOLUTELY NO WARRANTY. It is free software, and you are welcome to redistribute it under certain conditions (GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt). ABINIT is a project of the Universite Catholique de Louvain, Corning Inc. and other collaborators, see ~abinit/doc/developers/contributors.txt . Please read https://docs.abinit.org/theory/acknowledgments for suggested acknowledgments of the ABINIT effort. For more information, see https://www.abinit.org . .Starting date : Mon 24 Feb 2020. - ( at 16h44 ) - input file -> /home/gmatteo/git_repos/abinit/_abiref_gnu9.2_openmpi/tests/Test_suite/tutorial_tbase2_1/tbase2_1.in - output file -> tbase2_1.out - root for input files -> tbase2_1i - root for output files -> tbase2_1o DATASET 1 : space group P4/m m m (#123); Bravais tP (primitive tetrag.) ================================================================================ Values of the parameters that define the memory need for DATASET 1. intxc = 0 ionmov = 2 iscf = 7 lmnmax = 1 lnmax = 1 mgfft = 30 mpssoang = 1 mqgrid = 3001 natom = 2 nloc_mem = 1 nspden = 1 nspinor = 1 nsppol = 1 nsym = 16 n1xccc = 0 ntypat = 1 occopt = 1 xclevel = 1 - mband = 1 mffmem = 1 mkmem = 1 mpw = 752 nfft = 27000 nkpt = 1 ================================================================================ P This job should need less than 7.796 Mbytes of memory. Rough estimation (10% accuracy) of disk space for files : _ WF disk file : 0.013 Mbytes ; DEN or POT disk file : 0.208 Mbytes. ================================================================================ DATASET 2 : space group Pm -3 m (#221); Bravais cP (primitive cubic) ================================================================================ Values of the parameters that define the memory need for DATASET 2. intxc = 0 ionmov = 0 iscf = 7 lmnmax = 1 lnmax = 1 mgfft = 30 mpssoang = 1 mqgrid = 3001 natom = 1 nloc_mem = 1 nspden = 2 nspinor = 1 nsppol = 2 nsym = 48 n1xccc = 0 ntypat = 1 occopt = 2 xclevel = 1 - mband = 1 mffmem = 1 mkmem = 1 mpw = 752 nfft = 27000 nkpt = 1 ================================================================================ P This job should need less than 11.907 Mbytes of memory. Rough estimation (10% accuracy) of disk space for files : _ WF disk file : 0.025 Mbytes ; DEN or POT disk file : 0.414 Mbytes. ================================================================================ -------------------------------------------------------------------------------- ------------- Echo of variables that govern the present computation ------------ -------------------------------------------------------------------------------- - - outvars: echo of selected default values - iomode0 = 0 , fftalg0 =312 , wfoptalg0 = 0 - - outvars: echo of global parameters not present in the input file - max_nthreads = 0 - -outvars: echo values of preprocessed input variables -------- acell 1.0000000000E+01 1.0000000000E+01 1.0000000000E+01 Bohr amu 1.00794000E+00 bs_loband1 0 bs_loband2 0 0 diemac 2.00000000E+00 ecut 1.00000000E+01 Hartree - fftalg 312 ionmov1 2 ionmov2 0 istwfk 2 jdtset 1 2 kptopt 0 P mkmem 1 natom1 2 natom2 1 nband1 1 nband2 1 1 ndtset 2 ngfft 30 30 30 nkpt 1 nspden1 1 nspden2 2 nsppol1 1 nsppol2 2 nstep 10 nsym1 16 nsym2 48 ntime1 10 ntime2 1 ntypat 1 occ1 2.000000 occ2 1.000000 0.000000 occopt1 1 occopt2 2 optforces1 1 optforces2 2 spgroup1 123 spgroup2 221 spinat2 0.0000000000E+00 0.0000000000E+00 1.0000000000E+00 symafm1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 symafm2 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 symrel1 1 0 0 0 1 0 0 0 1 -1 0 0 0 -1 0 0 0 -1 -1 0 0 0 1 0 0 0 -1 1 0 0 0 -1 0 0 0 1 -1 0 0 0 -1 0 0 0 1 1 0 0 0 1 0 0 0 -1 1 0 0 0 -1 0 0 0 -1 -1 0 0 0 1 0 0 0 1 1 0 0 0 0 1 0 1 0 -1 0 0 0 0 -1 0 -1 0 -1 0 0 0 0 1 0 -1 0 1 0 0 0 0 -1 0 1 0 -1 0 0 0 0 -1 0 1 0 1 0 0 0 0 1 0 -1 0 1 0 0 0 0 -1 0 -1 0 -1 0 0 0 0 1 0 1 0 symrel2 1 0 0 0 1 0 0 0 1 -1 0 0 0 -1 0 0 0 -1 -1 0 0 0 1 0 0 0 -1 1 0 0 0 -1 0 0 0 1 -1 0 0 0 -1 0 0 0 1 1 0 0 0 1 0 0 0 -1 1 0 0 0 -1 0 0 0 -1 -1 0 0 0 1 0 0 0 1 0 1 0 1 0 0 0 0 1 0 -1 0 -1 0 0 0 0 -1 0 -1 0 1 0 0 0 0 -1 0 1 0 -1 0 0 0 0 1 0 -1 0 -1 0 0 0 0 1 0 1 0 1 0 0 0 0 -1 0 1 0 -1 0 0 0 0 -1 0 -1 0 1 0 0 0 0 1 0 0 1 1 0 0 0 1 0 0 0 -1 -1 0 0 0 -1 0 0 0 -1 1 0 0 0 -1 0 0 0 1 -1 0 0 0 1 0 0 0 -1 -1 0 0 0 1 0 0 0 1 1 0 0 0 -1 0 0 0 1 -1 0 0 0 -1 0 0 0 -1 1 0 0 0 1 0 1 0 0 0 0 1 0 1 0 -1 0 0 0 0 -1 0 -1 0 -1 0 0 0 0 1 0 -1 0 1 0 0 0 0 -1 0 1 0 -1 0 0 0 0 -1 0 1 0 1 0 0 0 0 1 0 -1 0 1 0 0 0 0 -1 0 -1 0 -1 0 0 0 0 1 0 1 0 0 1 0 0 0 1 1 0 0 0 -1 0 0 0 -1 -1 0 0 0 -1 0 0 0 1 -1 0 0 0 1 0 0 0 -1 1 0 0 0 -1 0 0 0 -1 1 0 0 0 1 0 0 0 1 -1 0 0 0 1 0 0 0 -1 -1 0 0 0 -1 0 0 0 1 1 0 0 0 0 1 0 1 0 1 0 0 0 0 -1 0 -1 0 -1 0 0 0 0 -1 0 1 0 -1 0 0 0 0 1 0 -1 0 1 0 0 0 0 -1 0 -1 0 1 0 0 0 0 1 0 1 0 -1 0 0 0 0 1 0 -1 0 -1 0 0 0 0 -1 0 1 0 1 0 0 tnons1 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 tnons2 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 toldfe1 0.00000000E+00 Hartree toldfe2 1.00000000E-06 Hartree toldff1 5.00000000E-05 toldff2 0.00000000E+00 tolmxf1 5.00000000E-04 tolmxf2 5.00000000E-05 typat1 1 1 typat2 1 xangst1 -3.7042404601E-01 0.0000000000E+00 0.0000000000E+00 3.7042404601E-01 0.0000000000E+00 0.0000000000E+00 xangst2 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 xcart1 -7.0000000000E-01 0.0000000000E+00 0.0000000000E+00 7.0000000000E-01 0.0000000000E+00 0.0000000000E+00 xcart2 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 xred1 -7.0000000000E-02 0.0000000000E+00 0.0000000000E+00 7.0000000000E-02 0.0000000000E+00 0.0000000000E+00 xred2 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 znucl 1.00000 ================================================================================ chkinp: Checking input parameters for consistency, jdtset= 1. chkinp: Checking input parameters for consistency, jdtset= 2. ================================================================================ == DATASET 1 ================================================================== - mpi_nproc: 1, omp_nthreads: -1 (-1 if OMP is not activated) --- !DatasetInfo iteration_state: {dtset: 1, } dimensions: {natom: 2, nkpt: 1, mband: 1, nsppol: 1, nspinor: 1, nspden: 1, mpw: 752, } cutoff_energies: {ecut: 10.0, pawecutdg: -1.0, } electrons: {nelect: 2.00000000E+00, charge: 0.00000000E+00, occopt: 1.00000000E+00, tsmear: 1.00000000E-02, } meta: {optdriver: 0, ionmov: 2, optcell: 0, iscf: 7, paral_kgb: 0, } ... Exchange-correlation functional for the present dataset will be: LDA: new Teter (4/93) with spin-polarized option - ixc=1 Citation for XC functional: S. Goedecker, M. Teter, J. Huetter, PRB 54, 1703 (1996) Real(R)+Recip(G) space primitive vectors, cartesian coordinates (Bohr,Bohr^-1): R(1)= 10.0000000 0.0000000 0.0000000 G(1)= 0.1000000 0.0000000 0.0000000 R(2)= 0.0000000 10.0000000 0.0000000 G(2)= 0.0000000 0.1000000 0.0000000 R(3)= 0.0000000 0.0000000 10.0000000 G(3)= 0.0000000 0.0000000 0.1000000 Unit cell volume ucvol= 1.0000000E+03 bohr^3 Angles (23,13,12)= 9.00000000E+01 9.00000000E+01 9.00000000E+01 degrees getcut: wavevector= 0.0000 0.0000 0.0000 ngfft= 30 30 30 ecut(hartree)= 10.000 => boxcut(ratio)= 2.10744 --- Pseudopotential description ------------------------------------------------ - pspini: atom type 1 psp file is /home/gmatteo/git_repos/abinit/tests/Psps_for_tests/01h.pspgth - pspatm: opening atomic psp file /home/gmatteo/git_repos/abinit/tests/Psps_for_tests/01h.pspgth - Goedecker-Teter-Hutter Wed May 8 14:27:44 EDT 1996 - 1.00000 1.00000 960508 znucl, zion, pspdat 2 1 0 0 2001 0.00000 pspcod,pspxc,lmax,lloc,mmax,r2well rloc= 0.2000000 cc1= -4.0663326; cc2= 0.6778322; cc3= 0.0000000; cc4= 0.0000000 rrs= 0.0000000; h1s= 0.0000000; h2s= 0.0000000 rrp= 0.0000000; h1p= 0.0000000 - Local part computed in reciprocal space. pspatm : COMMENT - the projectors are not normalized, so that the KB energies are not consistent with definition in PRB44, 8503 (1991). However, this does not influence the results obtained hereafter. pspatm : epsatm= -0.00480358 --- l ekb(1:nproj) --> pspatm: atomic psp has been read and splines computed -1.92143215E-02 ecore*ucvol(ha*bohr**3) -------------------------------------------------------------------------------- _setup2: Arith. and geom. avg. npw (full set) are 1503.000 1503.000 ================================================================================ === [ionmov= 2] Broyden-Fletcher-Goldfard-Shanno method (forces) ================================================================================ --- Iteration: ( 1/10) Internal Cycle: (1/1) -------------------------------------------------------------------------------- ---SELF-CONSISTENT-FIELD CONVERGENCE-------------------------------------------- --- !BeginCycle iteration_state: {dtset: 1, itime: 1, icycle: 1, } solver: {iscf: 7, nstep: 10, nline: 4, wfoptalg: 0, } tolerances: {toldff: 5.00E-05, } ... iter Etot(hartree) deltaE(h) residm vres2 diffor maxfor ETOT 1 -1.1013590048942 -1.101E+00 2.790E-06 8.389E+00 2.480E-02 2.480E-02 ETOT 2 -1.1036942492547 -2.335E-03 6.079E-10 2.843E-01 1.302E-02 3.782E-02 ETOT 3 -1.1037171066576 -2.286E-05 1.947E-08 1.545E-02 1.194E-03 3.662E-02 ETOT 4 -1.1037223545263 -5.248E-06 3.080E-08 2.641E-04 8.484E-04 3.747E-02 ETOT 5 -1.1037224211663 -6.664E-08 3.790E-10 8.035E-06 6.596E-05 3.740E-02 ETOT 6 -1.1037224213119 -1.456E-10 1.920E-13 3.970E-07 1.754E-06 3.741E-02 ETOT 7 -1.1037224213379 -2.606E-11 9.259E-14 9.779E-09 1.418E-06 3.740E-02 At SCF step 7, forces are converged : for the second time, max diff in force= 1.418E-06 < toldff= 5.000E-05 Cartesian components of stress tensor (hartree/bohr^3) sigma(1 1)= -1.64354976E-05 sigma(3 2)= 0.00000000E+00 sigma(2 2)= 3.60114991E-05 sigma(3 1)= 0.00000000E+00 sigma(3 3)= 3.60114991E-05 sigma(2 1)= 0.00000000E+00 --- !ResultsGS iteration_state: {dtset: 1, itime: 1, icycle: 1, } comment : Summary of ground state results lattice_vectors: - [ 10.0000000, 0.0000000, 0.0000000, ] - [ 0.0000000, 10.0000000, 0.0000000, ] - [ 0.0000000, 0.0000000, 10.0000000, ] lattice_lengths: [ 10.00000, 10.00000, 10.00000, ] lattice_angles: [ 90.000, 90.000, 90.000, ] # degrees, (23, 13, 12) lattice_volume: 1.0000000E+03 convergence: {deltae: -2.606E-11, res2: 9.779E-09, residm: 9.259E-14, diffor: 1.418E-06, } etotal : -1.10372242E+00 entropy : 0.00000000E+00 fermie : -3.65256341E-01 cartesian_stress_tensor: # hartree/bohr^3 - [ -1.64354976E-05, 0.00000000E+00, 0.00000000E+00, ] - [ 0.00000000E+00, 3.60114991E-05, 0.00000000E+00, ] - [ 0.00000000E+00, 0.00000000E+00, 3.60114991E-05, ] pressure_GPa: -5.4515E-01 xred : - [ -7.0000E-02, 0.0000E+00, 0.0000E+00, H] - [ 7.0000E-02, 0.0000E+00, 0.0000E+00, H] cartesian_forces: # hartree/bohr - [ -3.74042929E-02, -0.00000000E+00, -0.00000000E+00, ] - [ 3.74042929E-02, -0.00000000E+00, -0.00000000E+00, ] force_length_stats: {min: 3.74042929E-02, max: 3.74042929E-02, mean: 3.74042929E-02, } ... Integrated electronic density in atomic spheres: ------------------------------------------------ Atom Sphere_radius Integrated_density 1 2.00000 1.44122005 2 2.00000 1.44122005 ---OUTPUT----------------------------------------------------------------------- Cartesian coordinates (xcart) [bohr] -7.00000000000000E-01 0.00000000000000E+00 0.00000000000000E+00 7.00000000000000E-01 0.00000000000000E+00 0.00000000000000E+00 Reduced coordinates (xred) -7.00000000000000E-02 0.00000000000000E+00 0.00000000000000E+00 7.00000000000000E-02 0.00000000000000E+00 0.00000000000000E+00 Cartesian forces (fcart) [Ha/bohr]; max,rms= 3.74043E-02 2.15954E-02 (free atoms) -3.74042929459326E-02 -0.00000000000000E+00 -0.00000000000000E+00 3.74042929459326E-02 -0.00000000000000E+00 -0.00000000000000E+00 Reduced forces (fred) 3.74042929459326E-01 0.00000000000000E+00 0.00000000000000E+00 -3.74042929459326E-01 -0.00000000000000E+00 -0.00000000000000E+00 Total energy (etotal) [Ha]= -1.10372242133792E+00 --- Iteration: ( 2/10) Internal Cycle: (1/1) -------------------------------------------------------------------------------- ---SELF-CONSISTENT-FIELD CONVERGENCE-------------------------------------------- --- !BeginCycle iteration_state: {dtset: 1, itime: 2, icycle: 1, } solver: {iscf: 7, nstep: 10, nline: 4, wfoptalg: 0, } tolerances: {toldff: 5.00E-05, } ... iter Etot(hartree) deltaE(h) residm vres2 diffor maxfor ETOT 1 -1.1055395556802 -1.106E+00 2.497E-10 4.973E-02 2.480E-02 1.261E-02 ETOT 2 -1.1055522370018 -1.268E-05 3.023E-12 1.330E-03 2.024E-04 1.240E-02 ETOT 3 -1.1055525240950 -2.871E-07 2.475E-09 2.405E-04 2.101E-04 1.261E-02 ETOT 4 -1.1055525566531 -3.256E-08 1.770E-10 5.268E-06 5.187E-05 1.256E-02 ETOT 5 -1.1055525570209 -3.678E-10 3.130E-13 8.758E-08 2.693E-06 1.257E-02 ETOT 6 -1.1055525570220 -1.117E-12 9.726E-16 1.992E-09 9.891E-08 1.257E-02 At SCF step 6, forces are converged : for the second time, max diff in force= 9.891E-08 < toldff= 5.000E-05 Cartesian components of stress tensor (hartree/bohr^3) sigma(1 1)= 2.14754028E-05 sigma(3 2)= 0.00000000E+00 sigma(2 2)= 3.50015519E-05 sigma(3 1)= 0.00000000E+00 sigma(3 3)= 3.50015519E-05 sigma(2 1)= 0.00000000E+00 --- !ResultsGS iteration_state: {dtset: 1, itime: 2, icycle: 1, } comment : Summary of ground state results lattice_vectors: - [ 10.0000000, 0.0000000, 0.0000000, ] - [ 0.0000000, 10.0000000, 0.0000000, ] - [ 0.0000000, 0.0000000, 10.0000000, ] lattice_lengths: [ 10.00000, 10.00000, 10.00000, ] lattice_angles: [ 90.000, 90.000, 90.000, ] # degrees, (23, 13, 12) lattice_volume: 1.0000000E+03 convergence: {deltae: -1.117E-12, res2: 1.992E-09, residm: 9.726E-16, diffor: 9.891E-08, } etotal : -1.10555256E+00 entropy : 0.00000000E+00 fermie : -3.58935799E-01 cartesian_stress_tensor: # hartree/bohr^3 - [ 2.14754028E-05, 0.00000000E+00, 0.00000000E+00, ] - [ 0.00000000E+00, 3.50015519E-05, 0.00000000E+00, ] - [ 0.00000000E+00, 0.00000000E+00, 3.50015519E-05, ] pressure_GPa: -8.9713E-01 xred : - [ -7.3740E-02, 0.0000E+00, 0.0000E+00, H] - [ 7.3740E-02, 0.0000E+00, 0.0000E+00, H] cartesian_forces: # hartree/bohr - [ -1.25651505E-02, -0.00000000E+00, -0.00000000E+00, ] - [ 1.25651505E-02, -0.00000000E+00, -0.00000000E+00, ] force_length_stats: {min: 1.25651505E-02, max: 1.25651505E-02, mean: 1.25651505E-02, } ... Integrated electronic density in atomic spheres: ------------------------------------------------ Atom Sphere_radius Integrated_density 1 2.00000 1.40704092 2 2.00000 1.40704092 ---OUTPUT----------------------------------------------------------------------- Cartesian coordinates (xcart) [bohr] -7.37404292945932E-01 0.00000000000000E+00 0.00000000000000E+00 7.37404292945932E-01 0.00000000000000E+00 0.00000000000000E+00 Reduced coordinates (xred) -7.37404292945932E-02 0.00000000000000E+00 0.00000000000000E+00 7.37404292945932E-02 0.00000000000000E+00 0.00000000000000E+00 Cartesian forces (fcart) [Ha/bohr]; max,rms= 1.25652E-02 7.25449E-03 (free atoms) -1.25651504937955E-02 -0.00000000000000E+00 -0.00000000000000E+00 1.25651504937955E-02 -0.00000000000000E+00 -0.00000000000000E+00 Reduced forces (fred) 1.25651504937955E-01 0.00000000000000E+00 0.00000000000000E+00 -1.25651504937955E-01 -0.00000000000000E+00 -0.00000000000000E+00 Total energy (etotal) [Ha]= -1.10555255702201E+00 Difference of energy with previous step (new-old): Absolute (Ha)=-1.83014E-03 Relative =-1.65677E-03 --- Iteration: ( 3/10) Internal Cycle: (1/1) -------------------------------------------------------------------------------- ---SELF-CONSISTENT-FIELD CONVERGENCE-------------------------------------------- --- !BeginCycle iteration_state: {dtset: 1, itime: 3, icycle: 1, } solver: {iscf: 7, nstep: 10, nline: 4, wfoptalg: 0, } tolerances: {toldff: 5.00E-05, } ... iter Etot(hartree) deltaE(h) residm vres2 diffor maxfor ETOT 1 -1.1058239160568 -1.106E+00 6.373E-11 1.225E-02 1.039E-02 2.176E-03 ETOT 2 -1.1058269980950 -3.082E-06 7.141E-13 3.248E-04 9.796E-05 2.078E-03 ETOT 3 -1.1058270662905 -6.820E-08 5.827E-10 5.805E-05 9.995E-05 2.178E-03 ETOT 4 -1.1058270737818 -7.491E-09 4.047E-11 1.319E-06 2.398E-05 2.154E-03 ETOT 5 -1.1058270738704 -8.858E-11 7.625E-14 2.152E-08 1.309E-06 2.155E-03 At SCF step 5, forces are converged : for the second time, max diff in force= 1.309E-06 < toldff= 5.000E-05 Cartesian components of stress tensor (hartree/bohr^3) sigma(1 1)= 3.84829070E-05 sigma(3 2)= 0.00000000E+00 sigma(2 2)= 3.46517229E-05 sigma(3 1)= 0.00000000E+00 sigma(3 3)= 3.46517229E-05 sigma(2 1)= 0.00000000E+00 --- !ResultsGS iteration_state: {dtset: 1, itime: 3, icycle: 1, } comment : Summary of ground state results lattice_vectors: - [ 10.0000000, 0.0000000, 0.0000000, ] - [ 0.0000000, 10.0000000, 0.0000000, ] - [ 0.0000000, 0.0000000, 10.0000000, ] lattice_lengths: [ 10.00000, 10.00000, 10.00000, ] lattice_angles: [ 90.000, 90.000, 90.000, ] # degrees, (23, 13, 12) lattice_volume: 1.0000000E+03 convergence: {deltae: -8.858E-11, res2: 2.152E-08, residm: 7.625E-14, diffor: 1.309E-06, } etotal : -1.10582707E+00 entropy : 0.00000000E+00 fermie : -3.55897465E-01 cartesian_stress_tensor: # hartree/bohr^3 - [ 3.84829070E-05, 0.00000000E+00, 0.00000000E+00, ] - [ 0.00000000E+00, 3.46517229E-05, 0.00000000E+00, ] - [ 0.00000000E+00, 0.00000000E+00, 3.46517229E-05, ] pressure_GPa: -1.0571E+00 xred : - [ -7.5633E-02, 0.0000E+00, 0.0000E+00, H] - [ 7.5633E-02, 0.0000E+00, 0.0000E+00, H] cartesian_forces: # hartree/bohr - [ -2.15511232E-03, -0.00000000E+00, -0.00000000E+00, ] - [ 2.15511232E-03, -0.00000000E+00, -0.00000000E+00, ] force_length_stats: {min: 2.15511232E-03, max: 2.15511232E-03, mean: 2.15511232E-03, } ... Integrated electronic density in atomic spheres: ------------------------------------------------ Atom Sphere_radius Integrated_density 1 2.00000 1.39105245 2 2.00000 1.39105245 ---OUTPUT----------------------------------------------------------------------- Cartesian coordinates (xcart) [bohr] -7.56325661543310E-01 0.00000000000000E+00 0.00000000000000E+00 7.56325661543310E-01 0.00000000000000E+00 0.00000000000000E+00 Reduced coordinates (xred) -7.56325661543310E-02 0.00000000000000E+00 0.00000000000000E+00 7.56325661543310E-02 0.00000000000000E+00 0.00000000000000E+00 Cartesian forces (fcart) [Ha/bohr]; max,rms= 2.15511E-03 1.24425E-03 (free atoms) -2.15511231955776E-03 -0.00000000000000E+00 -0.00000000000000E+00 2.15511231955776E-03 -0.00000000000000E+00 -0.00000000000000E+00 Reduced forces (fred) 2.15511231955776E-02 0.00000000000000E+00 0.00000000000000E+00 -2.15511231955776E-02 -0.00000000000000E+00 -0.00000000000000E+00 Total energy (etotal) [Ha]= -1.10582707387041E+00 Difference of energy with previous step (new-old): Absolute (Ha)=-2.74517E-04 Relative =-2.48277E-04 --- Iteration: ( 4/10) Internal Cycle: (1/1) -------------------------------------------------------------------------------- ---SELF-CONSISTENT-FIELD CONVERGENCE-------------------------------------------- --- !BeginCycle iteration_state: {dtset: 1, itime: 4, icycle: 1, } solver: {iscf: 7, nstep: 10, nline: 4, wfoptalg: 0, } tolerances: {toldff: 5.00E-05, } ... iter Etot(hartree) deltaE(h) residm vres2 diffor maxfor ETOT 1 -1.1058359517763 -1.106E+00 2.639E-12 5.152E-04 1.995E-03 1.602E-04 ETOT 2 -1.1058360808934 -1.291E-07 2.916E-14 1.352E-05 2.273E-05 1.375E-04 ETOT 3 -1.1058360837580 -2.865E-09 2.432E-11 2.419E-06 2.042E-05 1.579E-04 At SCF step 3, forces are converged : for the second time, max diff in force= 2.042E-05 < toldff= 5.000E-05 Cartesian components of stress tensor (hartree/bohr^3) sigma(1 1)= 4.18328238E-05 sigma(3 2)= 0.00000000E+00 sigma(2 2)= 3.46056464E-05 sigma(3 1)= 0.00000000E+00 sigma(3 3)= 3.46056464E-05 sigma(2 1)= 0.00000000E+00 --- !ResultsGS iteration_state: {dtset: 1, itime: 4, icycle: 1, } comment : Summary of ground state results lattice_vectors: - [ 10.0000000, 0.0000000, 0.0000000, ] - [ 0.0000000, 10.0000000, 0.0000000, ] - [ 0.0000000, 0.0000000, 10.0000000, ] lattice_lengths: [ 10.00000, 10.00000, 10.00000, ] lattice_angles: [ 90.000, 90.000, 90.000, ] # degrees, (23, 13, 12) lattice_volume: 1.0000000E+03 convergence: {deltae: -2.865E-09, res2: 2.419E-06, residm: 2.432E-11, diffor: 2.042E-05, } etotal : -1.10583608E+00 entropy : 0.00000000E+00 fermie : -3.55256782E-01 cartesian_stress_tensor: # hartree/bohr^3 - [ 4.18328238E-05, 0.00000000E+00, 0.00000000E+00, ] - [ 0.00000000E+00, 3.46056464E-05, 0.00000000E+00, ] - [ 0.00000000E+00, 0.00000000E+00, 3.46056464E-05, ] pressure_GPa: -1.0890E+00 xred : - [ -7.6024E-02, 0.0000E+00, 0.0000E+00, H] - [ 7.6024E-02, 0.0000E+00, 0.0000E+00, H] cartesian_forces: # hartree/bohr - [ -1.57924482E-04, -0.00000000E+00, -0.00000000E+00, ] - [ 1.57924482E-04, -0.00000000E+00, -0.00000000E+00, ] force_length_stats: {min: 1.57924482E-04, max: 1.57924482E-04, mean: 1.57924482E-04, } ... Integrated electronic density in atomic spheres: ------------------------------------------------ Atom Sphere_radius Integrated_density 1 2.00000 1.38848194 2 2.00000 1.38848194 ---OUTPUT----------------------------------------------------------------------- Cartesian coordinates (xcart) [bohr] -7.60242810922072E-01 0.00000000000000E+00 0.00000000000000E+00 7.60242810922072E-01 0.00000000000000E+00 0.00000000000000E+00 Reduced coordinates (xred) -7.60242810922072E-02 0.00000000000000E+00 0.00000000000000E+00 7.60242810922072E-02 0.00000000000000E+00 0.00000000000000E+00 Cartesian forces (fcart) [Ha/bohr]; max,rms= 1.57924E-04 9.11777E-05 (free atoms) -1.57924482454065E-04 -0.00000000000000E+00 -0.00000000000000E+00 1.57924482454065E-04 -0.00000000000000E+00 -0.00000000000000E+00 Reduced forces (fred) 1.57924482454065E-03 0.00000000000000E+00 0.00000000000000E+00 -1.57924482454065E-03 -0.00000000000000E+00 -0.00000000000000E+00 Total energy (etotal) [Ha]= -1.10583608375795E+00 Difference of energy with previous step (new-old): Absolute (Ha)=-9.00989E-06 Relative =-8.14761E-06 At Broyd/MD step 4, gradients are converged : max grad (force/stress) = 1.5792E-04 < tolmxf= 5.0000E-04 ha/bohr (free atoms) ================================================================================ ----iterations are completed or convergence reached---- Mean square residual over all n,k,spin= 24.323E-12; max= 24.323E-12 reduced coordinates (array xred) for 2 atoms -0.076024281092 0.000000000000 0.000000000000 0.076024281092 0.000000000000 0.000000000000 rms dE/dt= 9.1178E-04; max dE/dt= 1.5792E-03; dE/dt below (all hartree) 1 0.001579244825 0.000000000000 0.000000000000 2 -0.001579244825 0.000000000000 0.000000000000 cartesian coordinates (angstrom) at end: 1 -0.40230316853436 0.00000000000000 0.00000000000000 2 0.40230316853436 0.00000000000000 0.00000000000000 cartesian forces (hartree/bohr) at end: 1 -0.00015792448245 -0.00000000000000 -0.00000000000000 2 0.00015792448245 -0.00000000000000 -0.00000000000000 frms,max,avg= 9.1177742E-05 1.5792448E-04 0.000E+00 0.000E+00 0.000E+00 h/b cartesian forces (eV/Angstrom) at end: 1 -0.00812080271635 -0.00000000000000 -0.00000000000000 2 0.00812080271635 -0.00000000000000 -0.00000000000000 frms,max,avg= 4.6885476E-03 8.1208027E-03 0.000E+00 0.000E+00 0.000E+00 e/A length scales= 10.000000000000 10.000000000000 10.000000000000 bohr = 5.291772085900 5.291772085900 5.291772085900 angstroms prteigrs : about to open file tbase2_1o_DS1_EIG Fermi (or HOMO) energy (hartree) = -0.35526 Average Vxc (hartree)= -0.07653 Eigenvalues (hartree) for nkpt= 1 k points: kpt# 1, nband= 1, wtk= 1.00000, kpt= 0.0000 0.0000 0.0000 (reduced coord) -0.35526 --- !EnergyTerms iteration_state : {dtset: 1, itime: 4, icycle: 1, } comment : Components of total free energy in Hartree kinetic : 9.50385455031023E-01 hartree : 6.78919676147409E-01 xc : -6.16674808436852E-01 Ewald energy : 9.52340313539042E-02 psp_core : -1.92143215271889E-05 local_psp : -2.21368122353191E+00 non_local_psp : 0.00000000000000E+00 total_energy : -1.10583608375795E+00 total_energy_eV : -3.00913301613768E+01 band_energy : -7.10513563405851E-01 ... rms coord change= 3.4781E-03 atom, delta coord (reduced): 1 -0.006024281092 0.000000000000 0.000000000000 2 0.006024281092 0.000000000000 0.000000000000 Cartesian components of stress tensor (hartree/bohr^3) sigma(1 1)= 4.18328238E-05 sigma(3 2)= 0.00000000E+00 sigma(2 2)= 3.46056464E-05 sigma(3 1)= 0.00000000E+00 sigma(3 3)= 3.46056464E-05 sigma(2 1)= 0.00000000E+00 -Cartesian components of stress tensor (GPa) [Pressure= -1.0890E+00 GPa] - sigma(1 1)= 1.23076396E+00 sigma(3 2)= 0.00000000E+00 - sigma(2 2)= 1.01813310E+00 sigma(3 1)= 0.00000000E+00 - sigma(3 3)= 1.01813310E+00 sigma(2 1)= 0.00000000E+00 ================================================================================ == DATASET 2 ================================================================== - mpi_nproc: 1, omp_nthreads: -1 (-1 if OMP is not activated) --- !DatasetInfo iteration_state: {dtset: 2, } dimensions: {natom: 1, nkpt: 1, mband: 1, nsppol: 2, nspinor: 1, nspden: 2, mpw: 752, } cutoff_energies: {ecut: 10.0, pawecutdg: -1.0, } electrons: {nelect: 1.00000000E+00, charge: 0.00000000E+00, occopt: 2.00000000E+00, tsmear: 1.00000000E-02, } meta: {optdriver: 0, ionmov: 0, optcell: 0, iscf: 7, paral_kgb: 0, } ... Exchange-correlation functional for the present dataset will be: LDA: new Teter (4/93) with spin-polarized option - ixc=1 Citation for XC functional: S. Goedecker, M. Teter, J. Huetter, PRB 54, 1703 (1996) Real(R)+Recip(G) space primitive vectors, cartesian coordinates (Bohr,Bohr^-1): R(1)= 10.0000000 0.0000000 0.0000000 G(1)= 0.1000000 0.0000000 0.0000000 R(2)= 0.0000000 10.0000000 0.0000000 G(2)= 0.0000000 0.1000000 0.0000000 R(3)= 0.0000000 0.0000000 10.0000000 G(3)= 0.0000000 0.0000000 0.1000000 Unit cell volume ucvol= 1.0000000E+03 bohr^3 Angles (23,13,12)= 9.00000000E+01 9.00000000E+01 9.00000000E+01 degrees getcut: wavevector= 0.0000 0.0000 0.0000 ngfft= 30 30 30 ecut(hartree)= 10.000 => boxcut(ratio)= 2.10744 -4.80358038E-03 ecore*ucvol(ha*bohr**3) -------------------------------------------------------------------------------- _setup2: Arith. and geom. avg. npw (full set) are 1503.000 1503.000 ================================================================================ --- !BeginCycle iteration_state: {dtset: 2, } solver: {iscf: 7, nstep: 10, nline: 4, wfoptalg: 0, } tolerances: {toldfe: 1.00E-06, } ... iter Etot(hartree) deltaE(h) residm vres2 ETOT 1 -0.46980757227245 -4.698E-01 5.815E-06 6.237E+00 ETOT 2 -0.47008829985619 -2.807E-04 2.519E-11 4.994E-01 ETOT 3 -0.47010514260040 -1.684E-05 1.440E-07 1.388E-02 ETOT 4 -0.47010529456329 -1.520E-07 2.587E-09 4.788E-04 ETOT 5 -0.47010531489229 -2.033E-08 1.467E-11 1.064E-05 At SCF step 5, etot is converged : for the second time, diff in etot= 2.033E-08 < toldfe= 1.000E-06 Cartesian components of stress tensor (hartree/bohr^3) sigma(1 1)= 1.80449486E-05 sigma(3 2)= 0.00000000E+00 sigma(2 2)= 1.80449486E-05 sigma(3 1)= 0.00000000E+00 sigma(3 3)= 1.80449486E-05 sigma(2 1)= 0.00000000E+00 --- !ResultsGS iteration_state: {dtset: 2, } comment : Summary of ground state results lattice_vectors: - [ 10.0000000, 0.0000000, 0.0000000, ] - [ 0.0000000, 10.0000000, 0.0000000, ] - [ 0.0000000, 0.0000000, 10.0000000, ] lattice_lengths: [ 10.00000, 10.00000, 10.00000, ] lattice_angles: [ 90.000, 90.000, 90.000, ] # degrees, (23, 13, 12) lattice_volume: 1.0000000E+03 convergence: {deltae: -2.033E-08, res2: 1.064E-05, residm: 1.467E-11, diffor: null, } etotal : -4.70105315E-01 entropy : 0.00000000E+00 fermie : -2.64144617E-01 cartesian_stress_tensor: # hartree/bohr^3 - [ 1.80449486E-05, 0.00000000E+00, 0.00000000E+00, ] - [ 0.00000000E+00, 1.80449486E-05, 0.00000000E+00, ] - [ 0.00000000E+00, 0.00000000E+00, 1.80449486E-05, ] pressure_GPa: -5.3090E-01 xred : - [ 0.0000E+00, 0.0000E+00, 0.0000E+00, H] cartesian_forces: # hartree/bohr - [ -0.00000000E+00, -0.00000000E+00, -0.00000000E+00, ] force_length_stats: {min: 0.00000000E+00, max: 0.00000000E+00, mean: 0.00000000E+00, } ... Integrated electronic and magnetization densities in atomic spheres: --------------------------------------------------------------------- Radius=ratsph(iatom), smearing ratsm= 0.0000. Diff(up-dn)=approximate z local magnetic moment. Atom Radius up_density dn_density Total(up+dn) Diff(up-dn) 1 2.00000 0.700509 0.000000 0.700509 0.700509 --------------------------------------------------------------------- Sum: 0.700509 0.000000 0.700509 0.700509 Total magnetization (from the atomic spheres): 0.700509 Total magnetization (exact up - dn): 1.000000 ================================================================================ ----iterations are completed or convergence reached---- Mean square residual over all n,k,spin= 14.158E-12; max= 14.666E-12 reduced coordinates (array xred) for 1 atoms 0.000000000000 0.000000000000 0.000000000000 rms dE/dt= 0.0000E+00; max dE/dt= 0.0000E+00; dE/dt below (all hartree) 1 0.000000000000 0.000000000000 0.000000000000 cartesian coordinates (angstrom) at end: 1 0.00000000000000 0.00000000000000 0.00000000000000 cartesian forces (hartree/bohr) at end: 1 -0.00000000000000 -0.00000000000000 -0.00000000000000 frms,max,avg= 0.0000000E+00 0.0000000E+00 0.000E+00 0.000E+00 0.000E+00 h/b cartesian forces (eV/Angstrom) at end: 1 -0.00000000000000 -0.00000000000000 -0.00000000000000 frms,max,avg= 0.0000000E+00 0.0000000E+00 0.000E+00 0.000E+00 0.000E+00 e/A length scales= 10.000000000000 10.000000000000 10.000000000000 bohr = 5.291772085900 5.291772085900 5.291772085900 angstroms prteigrs : about to open file tbase2_1o_DS2_EIG Fermi (or HOMO) energy (hartree) = -0.26414 Average Vxc (hartree)= -0.06898 Eigenvalues (hartree) for nkpt= 1 k points, SPIN UP: kpt# 1, nband= 1, wtk= 1.00000, kpt= 0.0000 0.0000 0.0000 (reduced coord) -0.26414 Eigenvalues (hartree) for nkpt= 1 k points, SPIN DOWN: kpt# 1, nband= 1, wtk= 1.00000, kpt= 0.0000 0.0000 0.0000 (reduced coord) -0.11112 --- !EnergyTerms iteration_state : {dtset: 2, } comment : Components of total free energy in Hartree kinetic : 4.06701917103141E-01 hartree : 1.47948764931797E-01 xc : -2.63312420910805E-01 Ewald energy : -1.41864873974033E-01 psp_core : -4.80358038179723E-06 local_psp : -6.19573898462006E-01 non_local_psp : 0.00000000000000E+00 total_energy : -4.70105314892287E-01 total_energy_eV : -1.27922161781602E+01 band_energy : -2.64144617221622E-01 ... Cartesian components of stress tensor (hartree/bohr^3) sigma(1 1)= 1.80449486E-05 sigma(3 2)= 0.00000000E+00 sigma(2 2)= 1.80449486E-05 sigma(3 1)= 0.00000000E+00 sigma(3 3)= 1.80449486E-05 sigma(2 1)= 0.00000000E+00 -Cartesian components of stress tensor (GPa) [Pressure= -5.3090E-01 GPa] - sigma(1 1)= 5.30900628E-01 sigma(3 2)= 0.00000000E+00 - sigma(2 2)= 5.30900628E-01 sigma(3 1)= 0.00000000E+00 - sigma(3 3)= 5.30900628E-01 sigma(2 1)= 0.00000000E+00 == END DATASET(S) ============================================================== ================================================================================ -outvars: echo values of variables after computation -------- acell 1.0000000000E+01 1.0000000000E+01 1.0000000000E+01 Bohr amu 1.00794000E+00 bs_loband1 0 bs_loband2 0 0 diemac 2.00000000E+00 ecut 1.00000000E+01 Hartree etotal1 -1.1058360838E+00 etotal2 -4.7010531489E-01 fcart1 -1.5792448245E-04 -0.0000000000E+00 -0.0000000000E+00 1.5792448245E-04 -0.0000000000E+00 -0.0000000000E+00 fcart2 -0.0000000000E+00 -0.0000000000E+00 -0.0000000000E+00 - fftalg 312 ionmov1 2 ionmov2 0 istwfk 2 jdtset 1 2 kptopt 0 P mkmem 1 natom1 2 natom2 1 nband1 1 nband2 1 1 ndtset 2 ngfft 30 30 30 nkpt 1 nspden1 1 nspden2 2 nsppol1 1 nsppol2 2 nstep 10 nsym1 16 nsym2 48 ntime1 10 ntime2 1 ntypat 1 occ1 2.000000 occ2 1.000000 0.000000 occopt1 1 occopt2 2 optforces1 1 optforces2 2 spgroup1 123 spgroup2 221 spinat2 0.0000000000E+00 0.0000000000E+00 1.0000000000E+00 strten1 4.1832823829E-05 3.4605646402E-05 3.4605646402E-05 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 strten2 1.8044948612E-05 1.8044948612E-05 1.8044948612E-05 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 symafm1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 symafm2 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 symrel1 1 0 0 0 1 0 0 0 1 -1 0 0 0 -1 0 0 0 -1 -1 0 0 0 1 0 0 0 -1 1 0 0 0 -1 0 0 0 1 -1 0 0 0 -1 0 0 0 1 1 0 0 0 1 0 0 0 -1 1 0 0 0 -1 0 0 0 -1 -1 0 0 0 1 0 0 0 1 1 0 0 0 0 1 0 1 0 -1 0 0 0 0 -1 0 -1 0 -1 0 0 0 0 1 0 -1 0 1 0 0 0 0 -1 0 1 0 -1 0 0 0 0 -1 0 1 0 1 0 0 0 0 1 0 -1 0 1 0 0 0 0 -1 0 -1 0 -1 0 0 0 0 1 0 1 0 symrel2 1 0 0 0 1 0 0 0 1 -1 0 0 0 -1 0 0 0 -1 -1 0 0 0 1 0 0 0 -1 1 0 0 0 -1 0 0 0 1 -1 0 0 0 -1 0 0 0 1 1 0 0 0 1 0 0 0 -1 1 0 0 0 -1 0 0 0 -1 -1 0 0 0 1 0 0 0 1 0 1 0 1 0 0 0 0 1 0 -1 0 -1 0 0 0 0 -1 0 -1 0 1 0 0 0 0 -1 0 1 0 -1 0 0 0 0 1 0 -1 0 -1 0 0 0 0 1 0 1 0 1 0 0 0 0 -1 0 1 0 -1 0 0 0 0 -1 0 -1 0 1 0 0 0 0 1 0 0 1 1 0 0 0 1 0 0 0 -1 -1 0 0 0 -1 0 0 0 -1 1 0 0 0 -1 0 0 0 1 -1 0 0 0 1 0 0 0 -1 -1 0 0 0 1 0 0 0 1 1 0 0 0 -1 0 0 0 1 -1 0 0 0 -1 0 0 0 -1 1 0 0 0 1 0 1 0 0 0 0 1 0 1 0 -1 0 0 0 0 -1 0 -1 0 -1 0 0 0 0 1 0 -1 0 1 0 0 0 0 -1 0 1 0 -1 0 0 0 0 -1 0 1 0 1 0 0 0 0 1 0 -1 0 1 0 0 0 0 -1 0 -1 0 -1 0 0 0 0 1 0 1 0 0 1 0 0 0 1 1 0 0 0 -1 0 0 0 -1 -1 0 0 0 -1 0 0 0 1 -1 0 0 0 1 0 0 0 -1 1 0 0 0 -1 0 0 0 -1 1 0 0 0 1 0 0 0 1 -1 0 0 0 1 0 0 0 -1 -1 0 0 0 -1 0 0 0 1 1 0 0 0 0 1 0 1 0 1 0 0 0 0 -1 0 -1 0 -1 0 0 0 0 -1 0 1 0 -1 0 0 0 0 1 0 -1 0 1 0 0 0 0 -1 0 -1 0 1 0 0 0 0 1 0 1 0 -1 0 0 0 0 1 0 -1 0 -1 0 0 0 0 -1 0 1 0 1 0 0 tnons1 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 tnons2 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 toldfe1 0.00000000E+00 Hartree toldfe2 1.00000000E-06 Hartree toldff1 5.00000000E-05 toldff2 0.00000000E+00 tolmxf1 5.00000000E-04 tolmxf2 5.00000000E-05 typat1 1 1 typat2 1 xangst1 -4.0230316853E-01 0.0000000000E+00 0.0000000000E+00 4.0230316853E-01 0.0000000000E+00 0.0000000000E+00 xangst2 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 xcart1 -7.6024281092E-01 0.0000000000E+00 0.0000000000E+00 7.6024281092E-01 0.0000000000E+00 0.0000000000E+00 xcart2 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 xred1 -7.6024281092E-02 0.0000000000E+00 0.0000000000E+00 7.6024281092E-02 0.0000000000E+00 0.0000000000E+00 xred2 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 znucl 1.00000 ================================================================================ - Timing analysis has been suppressed with timopt=0 ================================================================================ Suggested references for the acknowledgment of ABINIT usage. The users of ABINIT have little formal obligations with respect to the ABINIT group (those specified in the GNU General Public License, http://www.gnu.org/copyleft/gpl.txt). However, it is common practice in the scientific literature, to acknowledge the efforts of people that have made the research possible. In this spirit, please find below suggested citations of work written by ABINIT developers, corresponding to implementations inside of ABINIT that you have used in the present run. Note also that it will be of great value to readers of publications presenting these results, to read papers enabling them to understand the theoretical formalism and details of the ABINIT implementation. For information on why they are suggested, see also https://docs.abinit.org/theory/acknowledgments. - - [1] The Abinit project: Impact, environment and recent developments. - Computer Phys. Comm. 248, 107042 (2020). - X.Gonze, B. Amadon, G. Antonius, F.Arnardi, L.Baguet, J.-M.Beuken, - J.Bieder, F.Bottin, J.Bouchet, E.Bousquet, N.Brouwer, F.Bruneval, - G.Brunin, T.Cavignac, J.-B. Charraud, Wei Chen, M.Cote, S.Cottenier, - J.Denier, G.Geneste, Ph.Ghosez, M.Giantomassi, Y.Gillet, O.Gingras, - D.R.Hamann, G.Hautier, Xu He, N.Helbig, N.Holzwarth, Y.Jia, F.Jollet, - W.Lafargue-Dit-Hauret, K.Lejaeghere, M.A.L.Marques, A.Martin, C.Martins, - H.P.C. Miranda, F.Naccarato, K. Persson, G.Petretto, V.Planes, Y.Pouillon, - S.Prokhorenko, F.Ricci, G.-M.Rignanese, A.H.Romero, M.M.Schmitt, M.Torrent, - M.J.van Setten, B.Van Troeye, M.J.Verstraete, G.Zerah and J.W.Zwanzig - Comment: the fifth generic paper describing the ABINIT project. - Note that a version of this paper, that is not formatted for Computer Phys. Comm. - is available at https://www.abinit.org/sites/default/files/ABINIT20.pdf . - The licence allows the authors to put it on the Web. - DOI and bibtex: see https://docs.abinit.org/theory/bibliography/#gonze2020 - - [2] Recent developments in the ABINIT software package. - Computer Phys. Comm. 205, 106 (2016). - X.Gonze, F.Jollet, F.Abreu Araujo, D.Adams, B.Amadon, T.Applencourt, - C.Audouze, J.-M.Beuken, J.Bieder, A.Bokhanchuk, E.Bousquet, F.Bruneval - D.Caliste, M.Cote, F.Dahm, F.Da Pieve, M.Delaveau, M.Di Gennaro, - B.Dorado, C.Espejo, G.Geneste, L.Genovese, A.Gerossier, M.Giantomassi, - Y.Gillet, D.R.Hamann, L.He, G.Jomard, J.Laflamme Janssen, S.Le Roux, - A.Levitt, A.Lherbier, F.Liu, I.Lukacevic, A.Martin, C.Martins, - M.J.T.Oliveira, S.Ponce, Y.Pouillon, T.Rangel, G.-M.Rignanese, - A.H.Romero, B.Rousseau, O.Rubel, A.A.Shukri, M.Stankovski, M.Torrent, - M.J.Van Setten, B.Van Troeye, M.J.Verstraete, D.Waroquier, J.Wiktor, - B.Xu, A.Zhou, J.W.Zwanziger. - Comment: the fourth generic paper describing the ABINIT project. - Note that a version of this paper, that is not formatted for Computer Phys. Comm. - is available at https://www.abinit.org/sites/default/files/ABINIT16.pdf . - The licence allows the authors to put it on the Web. - DOI and bibtex: see https://docs.abinit.org/theory/bibliography/#gonze2016 - - [3] ABINIT: First-principles approach of materials and nanosystem properties. - Computer Phys. Comm. 180, 2582-2615 (2009). - X. Gonze, B. Amadon, P.-M. Anglade, J.-M. Beuken, F. Bottin, P. Boulanger, F. Bruneval, - D. Caliste, R. Caracas, M. Cote, T. Deutsch, L. Genovese, Ph. Ghosez, M. Giantomassi - S. Goedecker, D.R. Hamann, P. Hermet, F. Jollet, G. Jomard, S. Leroux, M. Mancini, S. Mazevet, - M.J.T. Oliveira, G. Onida, Y. Pouillon, T. Rangel, G.-M. Rignanese, D. Sangalli, R. Shaltaf, - M. Torrent, M.J. Verstraete, G. Zerah, J.W. Zwanziger - Comment: the third generic paper describing the ABINIT project. - Note that a version of this paper, that is not formatted for Computer Phys. Comm. - is available at https://www.abinit.org/sites/default/files/ABINIT_CPC_v10.pdf . - The licence allows the authors to put it on the Web. - DOI and bibtex: see https://docs.abinit.org/theory/bibliography/#gonze2009 - - And optionally: - - [4] A brief introduction to the ABINIT software package. - Z. Kristallogr. 220, 558-562 (2005). - X. Gonze, G.-M. Rignanese, M. Verstraete, J.-M. Beuken, Y. Pouillon, R. Caracas, F. Jollet, - M. Torrent, G. Zerah, M. Mikami, Ph. Ghosez, M. Veithen, J.-Y. Raty, V. Olevano, F. Bruneval, - L. Reining, R. Godby, G. Onida, D.R. Hamann, and D.C. Allan. - Comment: the second generic paper describing the ABINIT project. Note that this paper - should be cited especially if you are using the GW part of ABINIT, as several authors - of this part are not in the list of authors of the first or third paper. - The .pdf of the latter paper is available at https://www.abinit.org/sites/default/files/zfk_0505-06_558-562.pdf. - Note that it should not redistributed (Copyright by Oldenburg Wissenschaftverlag, - the licence allows the authors to put it on the Web). - DOI and bibtex: see https://docs.abinit.org/theory/bibliography/#gonze2005 - - Proc. 0 individual time (sec): cpu= 1.5 wall= 1.6 ================================================================================ Calculation completed. .Delivered 7 WARNINGs and 8 COMMENTs to log file. +Overall time at end (sec) : cpu= 1.5 wall= 1.6
You might use $ABI_TESTS/tutorial/Input/tbase2_x.files as files file (do not forget to modify it, like in tutorial 1, although it does not differ from tbase1_x.files.
../tbase2_1.in tbase2_x.out tbase2_xi tbase2_xo tbase2_x ../../../Psps_for_tests/01h.pspgth
Execute the code with:
abinit < tbase2_x.files > log 2> err &
The run should take less than one minute.
You should obtain the values:
etotal1 -1.1058360644E+00 etotal2 -4.7010531489E-01
and
xcart1 -7.6091015760E-01 0.0000000000E+00 0.0000000000E+00 7.6091015760E-01 0.0000000000E+00 0.0000000000E+00
These are similar to those determined in tutorial 1,
although they have been obtained in one run.
You can also check that the residual forces are lower than 5.0d-4
.
Convergence issues are discussed in section 7 of the abinit help file.
You should read it.
By the way, you have read many parts of the abinit help file!
You are missing the sections 2, 5, 7.
You are also missing the description of many input variables. We suggest that you finish reading entirely the abinit help file now, while the knowledge of the input variables will come in the long run.
2.1.b Many convergence parameters have already been identified. We will focus only on ecut and acell. This is because
-
the convergence of the SCF cycle and geometry determination are well under control thanks to toldfe, toldff and tolmxf (this might not be the case for other physical properties)
-
there is no k point convergence study to be done for an isolated system in a big box: no additional information is gained by adding a k-point beyond one
-
the boxcut value is automatically chosen larger than 2 by ABINIT, see the determination of the input variable ngfft by preprocessing
-
we are using ionmov = 2 for the determination of the geometry.
3 The convergence in ecut (II)¶
For the check of convergence with respect to ecut, you have the choice between doing different runs of the tbase2_1.in file with different values of ecut, or doing a double loop of datasets, as proposed in $ABI_TESTS/tutorial/Input/tbase2_2.in. The values of ecut have been chosen between 10 Ha and 35 Ha, by step of 5 Ha. If you want to make a double loop, you might benefit of reading again the double-loop section of the abinit_help file.
2.2.a You have likely seen a big increase of the CPU time needed to do the calculation. You should also look at the increase of the memory needed to do the calculation (go back to the beginning of the output file). The output data are as follows:
etotal11 -1.1058360644E+00 etotal12 -4.7010531489E-01 etotal21 -1.1218716100E+00 etotal22 -4.7529731401E-01 etotal31 -1.1291943792E+00 etotal32 -4.7773586424E-01 etotal41 -1.1326879404E+00 etotal42 -4.7899908214E-01 etotal51 -1.1346739190E+00 etotal52 -4.7972721653E-01 etotal61 -1.1359660026E+00 etotal62 -4.8022016459E-01 xcart11 -7.6091015760E-01 0.0000000000E+00 0.0000000000E+00 7.6091015760E-01 0.0000000000E+00 0.0000000000E+00 xcart12 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 xcart21 -7.5104912643E-01 0.0000000000E+00 0.0000000000E+00 7.5104912643E-01 0.0000000000E+00 0.0000000000E+00 xcart22 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 xcart31 -7.3977108831E-01 0.0000000000E+00 0.0000000000E+00 7.3977108831E-01 0.0000000000E+00 0.0000000000E+00 xcart32 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 xcart41 -7.3304273322E-01 0.0000000000E+00 0.0000000000E+00 7.3304273322E-01 0.0000000000E+00 0.0000000000E+00 xcart42 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 xcart51 -7.3001570260E-01 0.0000000000E+00 0.0000000000E+00 7.3001570260E-01 0.0000000000E+00 0.0000000000E+00 xcart52 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 xcart61 -7.2955902118E-01 0.0000000000E+00 0.0000000000E+00 7.2955902118E-01 0.0000000000E+00 0.0000000000E+00 xcart62 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00
The corresponding atomisation energies and interatomic distances are:
ecut (Ha) | atomisation energy (Ha) | interatomic distance (Bohr) |
---|---|---|
10 | .1656 | 1.522 |
15 | .1713 | 1.502 |
20 | .1737 | 1.480 |
25 | .1747 | 1.466 |
30 | .1753 | 1.460 |
35 | .1756 | 1.459 |
In order to obtain 0.2% relative accuracy on the bond length or atomisation energy, one should use a kinetic cut-off energy of 30 Ha. We will keep in mind this value for the final run.
4 The convergence in acell¶
The same technique as for ecut should be now used for the convergence in acell.
We will explore acell starting from 8 8 8
to 18 18 18
, by step of 2 2 2
.
We keep ecut 10 for this study. Indeed, it is a rather general rule that there is
little cross-influence between the convergence of ecut and the convergence of acell.
The file $ABI_TESTS/tutorial/Input/tbase2_3.in can be used as an example.
# H2 molecule in a big box # # This file to optimize the H2 bond length, compute the associated total # energy, then to compute the total energy of the isolated H atom. # Here, a double loop has been used. # ndtset 12 udtset 6 2 #Definition of the unit cell and ecut, #for which one will have to make a convergence study acell:? 8 8 8 acell+? 2 2 2 ecut 10 #First dataset : find the optimal bond length of H2, and associated total energy natom?1 2 # There are two atoms ionmov?1 2 # Use the modified Broyden algorithm ntime?1 10 # Maximum number of Broyden "timesteps" tolmxf?1 5.0d-4 # Stopping criterion for the geometry optimization : when # the residual forces are less than tolmxf, the Broyden # algorithm can stop xcart?1 -0.7 0.0 0.0 # The starting values of the 0.7 0.0 0.0 # atomic coordinates toldff?1 5.0d-5 # Will stop the SCF cycle when, twice in a row, # the difference between two consecutive evaluations of # forces differ by less than toldff (in Hartree/Bohr) nband?1 1 # Just one band #Second dataset : get the total energy of the isolated atom natom?2 1 # There is one atom nsppol?2 2 # Spin-polarized calculation occopt?2 2 # Allow occupation numbers to be set by hand nband?2 1 1 # Number of bands for spin up and spin down occ?2 1.0 0.0 # Occupation numbers for spin up state and spin down state. toldfe?2 1.0d-6 # Will stop the SCF cycles when, twice in a row, # the difference between two consecutive evaluations # of total energy differ by less than toldfe (in Hartree) # This value is way too large for most realistic studies of materials xcart?2 0.0 0.0 0.0 # The atom is located at the origin spinat?2 0.0 0.0 1.0 # Initialisation of spin #rprim 1 0 0 0 1 0 0 0 1 # This line, defining orthogonal primitive vectors, # is commented, because it is precisely the default value of rprim #Definition of the atom types ntypat 1 # There is only one type of atom znucl 1 # The keyword "znucl" refers to the atomic number of the # possible type(s) of atom. The pseudopotential(s) # mentioned in the "files" file must correspond # to the type(s) of atom. Here, the only type is Hydrogen. #Definition of the atoms typat 1 1 # For the first dataset, both numbers will be read, # while for the second dataset, only one number will be read #Definition of the k-point grid kptopt 0 # Enter the k points manually nkpt 1 # Only one k point is needed for isolated system, # taken by default to be 0.0 0.0 0.0 #Definition of the SCF procedure nstep 10 # Maximal number of SCF cycles #toldfe is no more defined, as toldff is used above... diemac 2.0 # Although this is not mandatory, it is worth to # precondition the SCF cycle. The model dielectric # function used as the standard preconditioner # is described in the "dielng" input variable section. # Here, we follow the prescriptions for molecules # in a big box #%%<BEGIN TEST_INFO> #%% [setup] #%% executable = abinit #%% [files] #%% files_to_test = #%% tbase2_3.out, tolnlines= 0, tolabs= 0.000e+00, tolrel= 0.000e+00, fld_options=-medium #%% psp_files = 01h.pspgth #%% [paral_info] #%% max_nprocs = 1 #%% [extra_info] #%% authors = Unknown #%% keywords = #%% description = #%% H2 molecule in a big box #%% This file to optimize the H2 bond length, compute the associated total #%% energy, then to compute the total energy of the isolated H atom. #%% Here, a double loop has been used. #%%<END TEST_INFO>
The output results in $ABI_TESTS/tutorial/Refs/tbase2_3.out are as follows:
etotal11 -1.1188124709E+00 etotal12 -4.8074164402E-01 etotal21 -1.1058360838E+00 etotal22 -4.7010531489E-01 etotal31 -1.1039109527E+00 etotal32 -4.6767804802E-01 etotal41 -1.1039012868E+00 etotal42 -4.6743724199E-01 etotal51 -1.1041439411E+00 etotal52 -4.6735895176E-01 etotal61 -1.1042058281E+00 etotal62 -4.6736729718E-01 xcart11 -7.8330751426E-01 0.0000000000E+00 0.0000000000E+00 7.8330751426E-01 0.0000000000E+00 0.0000000000E+00 xcart12 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 xcart21 -7.6024281092E-01 0.0000000000E+00 0.0000000000E+00 7.6024281092E-01 0.0000000000E+00 0.0000000000E+00 xcart22 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 xcart31 -7.5428234893E-01 0.0000000000E+00 0.0000000000E+00 7.5428234893E-01 0.0000000000E+00 0.0000000000E+00 xcart32 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 xcart41 -7.5446921004E-01 0.0000000000E+00 0.0000000000E+00 7.5446921004E-01 0.0000000000E+00 0.0000000000E+00 xcart42 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 xcart51 -7.5384974520E-01 0.0000000000E+00 0.0000000000E+00 7.5384974520E-01 0.0000000000E+00 0.0000000000E+00 xcart52 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 xcart61 -7.5373336127E-01 0.0000000000E+00 0.0000000000E+00 7.5373336127E-01 0.0000000000E+00 0.0000000000E+00 xcart62 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00
The corresponding atomisation energies and interatomic distances are:
acell (Bohr) | atomisation energy (Ha) | interatomic distance (Bohr) |
---|---|---|
8 | .1574 | 1.568 |
10 | .1656 | 1.522 |
12 | .1686 | 1.509 |
14 | .1691 | 1.510 |
16 | .1694 | 1.508 |
18 | .1695 | 1.508 |
In order to reach 0.2% convergence on the interatomic distance, one needs acell 12 12 12
.
The atomisation energy needs acell 14 14 14
to be converged at that level.
At 12 12 12
, the difference is .0009 Ha = 0.024 eV, which is sufficiently small for practical purposes.
We will use acell 12 12 12
for the final run.
For most solids the size of the unit cell will be smaller than that. We are treating a lot of vacuum in this supercell! So, the H_2 study, with this pseudopotential, turns out to be not really easy. Of course, the number of states to be treated is minimal! This allows to have reasonable CPU time still.
5 The final calculation in Local (Spin) Density Approximation¶
We now use the correct values of both ecut and acell.
Well, you should modify the tbase2_3.in file to make a calculation with acell 12 12 12
and ecut 30
.
You can still use the double loop feature with udtset 1 2
(which reduces to a single loop), to minimize the modifications to the file.
The file $ABI_TESTS/tutorial/Input/tbase2_4.in can be taken as an example of input file:
# H2 molecule in a big box # # This file to optimize the H2 bond length, compute the associated total # energy, then to compute the total energy of the isolated H atom. # Here, the ecut and acell are fixed : the double loop reduces # effectively to a single loop. ndtset 2 udtset 1 2 #Definition of the unit cell and ecut, #for which one will have to make a convergence study acell:? 12 12 12 acell+? 2 2 2 ecut 30 #First dataset : find the optimal bond length of H2, and associated total energy natom?1 2 # There are two atoms ionmov?1 2 # Use the modified Broyden algorithm ntime?1 10 # Maximum number of Broyden "timesteps" tolmxf?1 5.0d-4 # Stopping criterion for the geometry optimization : when # the residual forces are less than tolmxf, the Broyden # algorithm can stop xcart?1 -0.7 0.0 0.0 # The starting values of the 0.7 0.0 0.0 # atomic coordinates toldff?1 5.0d-5 # Will stop the SCF cycle when, twice in a row, # the difference between two consecutive evaluations of # forces differ by less than toldff (in Hartree/Bohr) nband?1 1 # Just one band #Second dataset : get the total energy of the isolated atom natom?2 1 # There is one atom nsppol?2 2 # Spin-polarized calculation occopt?2 2 # Allow occupation numbers to be set by hand nband?2 1 1 # Number of bands for spin up and spin down occ?2 1.0 0.0 # Occupation numbers for spin up state and spin down state. toldfe?2 1.0d-6 # Will stop the SCF cycles when, twice in a row, # the difference between two consecutive evaluations # of total energy differ by less than toldfe (in Hartree) # This value is way too large for most realistic studies of materials xcart?2 0.0 0.0 0.0 # The atom is located at the origin spinat?2 0.0 0.0 1.0 # Initialisation of spin #rprim 1 0 0 0 1 0 0 0 1 # This line, defining orthogonal primitive vectors, # is commented, because it is precisely the default value of rprim #Definition of the atom types ntypat 1 # There is only one type of atom znucl 1 # The keyword "znucl" refers to the atomic number of the # possible type(s) of atom. The pseudopotential(s) # mentioned in the "files" file must correspond # to the type(s) of atom. Here, the only type is Hydrogen. #Definition of the atoms typat 1 1 # For the first dataset, both numbers will be read, # while for the second dataset, only one number will be read #Definition of the k-point grid kptopt 0 # Enter the k points manually nkpt 1 # Only one k point is needed for isolated system, # taken by default to be 0.0 0.0 0.0 #Definition of the SCF procedure nstep 10 # Maximal number of SCF cycles #toldfe is no more defined, as toldff is used above... diemac 1.0 # Although this is not mandatory, it is worth to # precondition the SCF cycle. The model dielectric # function used as the standard preconditioner # is described in the "dielng" input variable section. # Here, we follow the prescriptions for molecules # in a big box #%%<BEGIN TEST_INFO> #%% [setup] #%% executable = abinit #%% [files] #%% files_to_test = #%% tbase2_4.out, tolnlines= 0, tolabs= 0.000e+00, tolrel= 0.000e+00, fld_options=-medium #%% psp_files = 01h.pspgth #%% [paral_info] #%% max_nprocs = 1 #%% [extra_info] #%% authors = Unknown #%% keywords = #%% description = #%% H2 molecule in a big box #%% This file to optimize the H2 bond length, compute the associated total #%% energy, then to compute the total energy of the isolated H atom. #%% Here, the ecut and acell are fixed : the double loop reduces #%% effectively to a single loop. #%%<END TEST_INFO>
while $ABI_TESTS/tutorial/Refs/tbase2_4.out is as an example of output file:
.Version 9.0.0 of ABINIT .(MPI version, prepared for a x86_64_linux_gnu9.2 computer) .Copyright (C) 1998-2020 ABINIT group . ABINIT comes with ABSOLUTELY NO WARRANTY. It is free software, and you are welcome to redistribute it under certain conditions (GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt). ABINIT is a project of the Universite Catholique de Louvain, Corning Inc. and other collaborators, see ~abinit/doc/developers/contributors.txt . Please read https://docs.abinit.org/theory/acknowledgments for suggested acknowledgments of the ABINIT effort. For more information, see https://www.abinit.org . .Starting date : Mon 24 Feb 2020. - ( at 16h44 ) - input file -> /home/gmatteo/git_repos/abinit/_abiref_gnu9.2_openmpi/tests/Test_suite/tutorial_tbase2_4/tbase2_4.in - output file -> tbase2_4.out - root for input files -> tbase2_4i - root for output files -> tbase2_4o DATASET 11 : space group P4/m m m (#123); Bravais tP (primitive tetrag.) ================================================================================ Values of the parameters that define the memory need for DATASET 11. intxc = 0 ionmov = 2 iscf = 7 lmnmax = 1 lnmax = 1 mgfft = 60 mpssoang = 1 mqgrid = 3001 natom = 2 nloc_mem = 1 nspden = 1 nspinor = 1 nsppol = 1 nsym = 16 n1xccc = 0 ntypat = 1 occopt = 1 xclevel = 1 - mband = 1 mffmem = 1 mkmem = 1 mpw = 6759 nfft = 216000 nkpt = 1 ================================================================================ P This job should need less than 59.627 Mbytes of memory. Rough estimation (10% accuracy) of disk space for files : _ WF disk file : 0.105 Mbytes ; DEN or POT disk file : 1.650 Mbytes. ================================================================================ DATASET 12 : space group Pm -3 m (#221); Bravais cP (primitive cubic) ================================================================================ Values of the parameters that define the memory need for DATASET 12. intxc = 0 ionmov = 0 iscf = 7 lmnmax = 1 lnmax = 1 mgfft = 60 mpssoang = 1 mqgrid = 3001 natom = 1 nloc_mem = 1 nspden = 2 nspinor = 1 nsppol = 2 nsym = 48 n1xccc = 0 ntypat = 1 occopt = 2 xclevel = 1 - mband = 1 mffmem = 1 mkmem = 1 mpw = 6759 nfft = 216000 nkpt = 1 ================================================================================ P This job should need less than 92.575 Mbytes of memory. Rough estimation (10% accuracy) of disk space for files : _ WF disk file : 0.208 Mbytes ; DEN or POT disk file : 3.298 Mbytes. ================================================================================ -------------------------------------------------------------------------------- ------------- Echo of variables that govern the present computation ------------ -------------------------------------------------------------------------------- - - outvars: echo of selected default values - iomode0 = 0 , fftalg0 =312 , wfoptalg0 = 0 - - outvars: echo of global parameters not present in the input file - max_nthreads = 0 - -outvars: echo values of preprocessed input variables -------- acell 1.2000000000E+01 1.2000000000E+01 1.2000000000E+01 Bohr amu 1.00794000E+00 bs_loband11 0 bs_loband12 0 0 diemac 1.00000000E+00 ecut 3.00000000E+01 Hartree - fftalg 312 ionmov11 2 ionmov12 0 istwfk 2 jdtset 11 12 kptopt 0 P mkmem 1 natom11 2 natom12 1 nband11 1 nband12 1 1 ndtset 2 ngfft 60 60 60 nkpt 1 nspden11 1 nspden12 2 nsppol11 1 nsppol12 2 nstep 10 nsym11 16 nsym12 48 ntime11 10 ntime12 1 ntypat 1 occ11 2.000000 occ12 1.000000 0.000000 occopt11 1 occopt12 2 optforces11 1 optforces12 2 spgroup11 123 spgroup12 221 spinat12 0.0000000000E+00 0.0000000000E+00 1.0000000000E+00 symafm11 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 symafm12 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 symrel11 1 0 0 0 1 0 0 0 1 -1 0 0 0 -1 0 0 0 -1 -1 0 0 0 1 0 0 0 -1 1 0 0 0 -1 0 0 0 1 -1 0 0 0 -1 0 0 0 1 1 0 0 0 1 0 0 0 -1 1 0 0 0 -1 0 0 0 -1 -1 0 0 0 1 0 0 0 1 1 0 0 0 0 1 0 1 0 -1 0 0 0 0 -1 0 -1 0 -1 0 0 0 0 1 0 -1 0 1 0 0 0 0 -1 0 1 0 -1 0 0 0 0 -1 0 1 0 1 0 0 0 0 1 0 -1 0 1 0 0 0 0 -1 0 -1 0 -1 0 0 0 0 1 0 1 0 symrel12 1 0 0 0 1 0 0 0 1 -1 0 0 0 -1 0 0 0 -1 -1 0 0 0 1 0 0 0 -1 1 0 0 0 -1 0 0 0 1 -1 0 0 0 -1 0 0 0 1 1 0 0 0 1 0 0 0 -1 1 0 0 0 -1 0 0 0 -1 -1 0 0 0 1 0 0 0 1 0 1 0 1 0 0 0 0 1 0 -1 0 -1 0 0 0 0 -1 0 -1 0 1 0 0 0 0 -1 0 1 0 -1 0 0 0 0 1 0 -1 0 -1 0 0 0 0 1 0 1 0 1 0 0 0 0 -1 0 1 0 -1 0 0 0 0 -1 0 -1 0 1 0 0 0 0 1 0 0 1 1 0 0 0 1 0 0 0 -1 -1 0 0 0 -1 0 0 0 -1 1 0 0 0 -1 0 0 0 1 -1 0 0 0 1 0 0 0 -1 -1 0 0 0 1 0 0 0 1 1 0 0 0 -1 0 0 0 1 -1 0 0 0 -1 0 0 0 -1 1 0 0 0 1 0 1 0 0 0 0 1 0 1 0 -1 0 0 0 0 -1 0 -1 0 -1 0 0 0 0 1 0 -1 0 1 0 0 0 0 -1 0 1 0 -1 0 0 0 0 -1 0 1 0 1 0 0 0 0 1 0 -1 0 1 0 0 0 0 -1 0 -1 0 -1 0 0 0 0 1 0 1 0 0 1 0 0 0 1 1 0 0 0 -1 0 0 0 -1 -1 0 0 0 -1 0 0 0 1 -1 0 0 0 1 0 0 0 -1 1 0 0 0 -1 0 0 0 -1 1 0 0 0 1 0 0 0 1 -1 0 0 0 1 0 0 0 -1 -1 0 0 0 -1 0 0 0 1 1 0 0 0 0 1 0 1 0 1 0 0 0 0 -1 0 -1 0 -1 0 0 0 0 -1 0 1 0 -1 0 0 0 0 1 0 -1 0 1 0 0 0 0 -1 0 -1 0 1 0 0 0 0 1 0 1 0 -1 0 0 0 0 1 0 -1 0 -1 0 0 0 0 -1 0 1 0 1 0 0 tnons11 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 tnons12 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 toldfe11 0.00000000E+00 Hartree toldfe12 1.00000000E-06 Hartree toldff11 5.00000000E-05 toldff12 0.00000000E+00 tolmxf11 5.00000000E-04 tolmxf12 5.00000000E-05 typat11 1 1 typat12 1 xangst11 -3.7042404601E-01 0.0000000000E+00 0.0000000000E+00 3.7042404601E-01 0.0000000000E+00 0.0000000000E+00 xangst12 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 xcart11 -7.0000000000E-01 0.0000000000E+00 0.0000000000E+00 7.0000000000E-01 0.0000000000E+00 0.0000000000E+00 xcart12 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 xred11 -5.8333333333E-02 0.0000000000E+00 0.0000000000E+00 5.8333333333E-02 0.0000000000E+00 0.0000000000E+00 xred12 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 znucl 1.00000 ================================================================================ chkinp: Checking input parameters for consistency, jdtset= 11. chkinp: Checking input parameters for consistency, jdtset= 12. ================================================================================ == DATASET 11 ================================================================== - mpi_nproc: 1, omp_nthreads: -1 (-1 if OMP is not activated) --- !DatasetInfo iteration_state: {dtset: 11, } dimensions: {natom: 2, nkpt: 1, mband: 1, nsppol: 1, nspinor: 1, nspden: 1, mpw: 6759, } cutoff_energies: {ecut: 30.0, pawecutdg: -1.0, } electrons: {nelect: 2.00000000E+00, charge: 0.00000000E+00, occopt: 1.00000000E+00, tsmear: 1.00000000E-02, } meta: {optdriver: 0, ionmov: 2, optcell: 0, iscf: 7, paral_kgb: 0, } ... Exchange-correlation functional for the present dataset will be: LDA: new Teter (4/93) with spin-polarized option - ixc=1 Citation for XC functional: S. Goedecker, M. Teter, J. Huetter, PRB 54, 1703 (1996) Real(R)+Recip(G) space primitive vectors, cartesian coordinates (Bohr,Bohr^-1): R(1)= 12.0000000 0.0000000 0.0000000 G(1)= 0.0833333 0.0000000 0.0000000 R(2)= 0.0000000 12.0000000 0.0000000 G(2)= 0.0000000 0.0833333 0.0000000 R(3)= 0.0000000 0.0000000 12.0000000 G(3)= 0.0000000 0.0000000 0.0833333 Unit cell volume ucvol= 1.7280000E+03 bohr^3 Angles (23,13,12)= 9.00000000E+01 9.00000000E+01 9.00000000E+01 degrees getcut: wavevector= 0.0000 0.0000 0.0000 ngfft= 60 60 60 ecut(hartree)= 30.000 => boxcut(ratio)= 2.02789 --- Pseudopotential description ------------------------------------------------ - pspini: atom type 1 psp file is /home/gmatteo/git_repos/abinit/tests/Psps_for_tests/01h.pspgth - pspatm: opening atomic psp file /home/gmatteo/git_repos/abinit/tests/Psps_for_tests/01h.pspgth - Goedecker-Teter-Hutter Wed May 8 14:27:44 EDT 1996 - 1.00000 1.00000 960508 znucl, zion, pspdat 2 1 0 0 2001 0.00000 pspcod,pspxc,lmax,lloc,mmax,r2well rloc= 0.2000000 cc1= -4.0663326; cc2= 0.6778322; cc3= 0.0000000; cc4= 0.0000000 rrs= 0.0000000; h1s= 0.0000000; h2s= 0.0000000 rrp= 0.0000000; h1p= 0.0000000 - Local part computed in reciprocal space. pspatm : COMMENT - the projectors are not normalized, so that the KB energies are not consistent with definition in PRB44, 8503 (1991). However, this does not influence the results obtained hereafter. pspatm : epsatm= -0.00480358 --- l ekb(1:nproj) --> pspatm: atomic psp has been read and splines computed -1.92143215E-02 ecore*ucvol(ha*bohr**3) -------------------------------------------------------------------------------- _setup2: Arith. and geom. avg. npw (full set) are 13517.000 13517.000 ================================================================================ === [ionmov= 2] Broyden-Fletcher-Goldfard-Shanno method (forces) ================================================================================ --- Iteration: ( 1/10) Internal Cycle: (1/1) -------------------------------------------------------------------------------- ---SELF-CONSISTENT-FIELD CONVERGENCE-------------------------------------------- --- !BeginCycle iteration_state: {dtset: 11, itime: 1, icycle: 1, } solver: {iscf: 7, nstep: 10, nline: 4, wfoptalg: 0, } tolerances: {toldff: 5.00E-05, } ... iter Etot(hartree) deltaE(h) residm vres2 diffor maxfor ETOT 1 -1.1294675706715 -1.129E+00 4.495E-05 4.802E+01 7.290E-03 7.290E-03 ETOT 2 -1.1320826829733 -2.615E-03 1.618E-09 3.564E+00 1.804E-02 2.533E-02 ETOT 3 -1.1324276037537 -3.449E-04 4.497E-07 1.037E-01 7.182E-03 1.815E-02 ETOT 4 -1.1324407193870 -1.312E-05 1.024E-07 8.543E-05 1.099E-03 1.925E-02 ETOT 5 -1.1324407238201 -4.433E-09 6.820E-12 1.579E-05 1.186E-05 1.926E-02 ETOT 6 -1.1324407238088 1.131E-11 2.743E-14 8.002E-07 1.745E-06 1.926E-02 At SCF step 6, forces are converged : for the second time, max diff in force= 1.745E-06 < toldff= 5.000E-05 Cartesian components of stress tensor (hartree/bohr^3) sigma(1 1)= -1.18576493E-05 sigma(3 2)= 0.00000000E+00 sigma(2 2)= 4.87250422E-06 sigma(3 1)= 0.00000000E+00 sigma(3 3)= 4.87250422E-06 sigma(2 1)= 0.00000000E+00 --- !ResultsGS iteration_state: {dtset: 11, itime: 1, icycle: 1, } comment : Summary of ground state results lattice_vectors: - [ 12.0000000, 0.0000000, 0.0000000, ] - [ 0.0000000, 12.0000000, 0.0000000, ] - [ 0.0000000, 0.0000000, 12.0000000, ] lattice_lengths: [ 12.00000, 12.00000, 12.00000, ] lattice_angles: [ 90.000, 90.000, 90.000, ] # degrees, (23, 13, 12) lattice_volume: 1.7280000E+03 convergence: {deltae: 1.131E-11, res2: 8.002E-07, residm: 2.743E-14, diffor: 1.745E-06, } etotal : -1.13244072E+00 entropy : 0.00000000E+00 fermie : -3.71601829E-01 cartesian_stress_tensor: # hartree/bohr^3 - [ -1.18576493E-05, 0.00000000E+00, 0.00000000E+00, ] - [ 0.00000000E+00, 4.87250422E-06, 0.00000000E+00, ] - [ 0.00000000E+00, 0.00000000E+00, 4.87250422E-06, ] pressure_GPa: 2.0719E-02 xred : - [ -5.8333E-02, 0.0000E+00, 0.0000E+00, H] - [ 5.8333E-02, 0.0000E+00, 0.0000E+00, H] cartesian_forces: # hartree/bohr - [ -1.92592034E-02, -0.00000000E+00, -0.00000000E+00, ] - [ 1.92592034E-02, -0.00000000E+00, -0.00000000E+00, ] force_length_stats: {min: 1.92592034E-02, max: 1.92592034E-02, mean: 1.92592034E-02, } ... Integrated electronic density in atomic spheres: ------------------------------------------------ Atom Sphere_radius Integrated_density 1 2.00000 1.47503350 2 2.00000 1.47503350 ---OUTPUT----------------------------------------------------------------------- Cartesian coordinates (xcart) [bohr] -7.00000000000000E-01 0.00000000000000E+00 0.00000000000000E+00 7.00000000000000E-01 0.00000000000000E+00 0.00000000000000E+00 Reduced coordinates (xred) -5.83333333333333E-02 0.00000000000000E+00 0.00000000000000E+00 5.83333333333333E-02 0.00000000000000E+00 0.00000000000000E+00 Cartesian forces (fcart) [Ha/bohr]; max,rms= 1.92592E-02 1.11193E-02 (free atoms) -1.92592033774355E-02 -0.00000000000000E+00 -0.00000000000000E+00 1.92592033774355E-02 -0.00000000000000E+00 -0.00000000000000E+00 Reduced forces (fred) 2.31110440529226E-01 0.00000000000000E+00 0.00000000000000E+00 -2.31110440529226E-01 -0.00000000000000E+00 -0.00000000000000E+00 Total energy (etotal) [Ha]= -1.13244072380879E+00 --- Iteration: ( 2/10) Internal Cycle: (1/1) -------------------------------------------------------------------------------- ---SELF-CONSISTENT-FIELD CONVERGENCE-------------------------------------------- --- !BeginCycle iteration_state: {dtset: 11, itime: 2, icycle: 1, } solver: {iscf: 7, nstep: 10, nline: 4, wfoptalg: 0, } tolerances: {toldff: 5.00E-05, } ... iter Etot(hartree) deltaE(h) residm vres2 diffor maxfor ETOT 1 -1.1328974456390 -1.133E+00 1.653E-10 9.284E-02 1.517E-02 4.088E-03 ETOT 2 -1.1329006332982 -3.188E-06 1.791E-12 3.456E-03 6.587E-04 4.746E-03 ETOT 3 -1.1329010352884 -4.020E-07 4.901E-10 1.175E-04 2.481E-04 4.995E-03 ETOT 4 -1.1329010475859 -1.230E-08 9.256E-11 1.288E-07 2.965E-05 4.965E-03 ETOT 5 -1.1329010475893 -3.411E-12 3.093E-15 1.698E-08 6.089E-07 4.964E-03 At SCF step 5, forces are converged : for the second time, max diff in force= 6.089E-07 < toldff= 5.000E-05 Cartesian components of stress tensor (hartree/bohr^3) sigma(1 1)= -4.02597615E-07 sigma(3 2)= 0.00000000E+00 sigma(2 2)= 4.75261789E-06 sigma(3 1)= 0.00000000E+00 sigma(3 3)= 4.75261789E-06 sigma(2 1)= 0.00000000E+00 --- !ResultsGS iteration_state: {dtset: 11, itime: 2, icycle: 1, } comment : Summary of ground state results lattice_vectors: - [ 12.0000000, 0.0000000, 0.0000000, ] - [ 0.0000000, 12.0000000, 0.0000000, ] - [ 0.0000000, 0.0000000, 12.0000000, ] lattice_lengths: [ 12.00000, 12.00000, 12.00000, ] lattice_angles: [ 90.000, 90.000, 90.000, ] # degrees, (23, 13, 12) lattice_volume: 1.7280000E+03 convergence: {deltae: -3.411E-12, res2: 1.698E-08, residm: 3.093E-15, diffor: 6.089E-07, } etotal : -1.13290105E+00 entropy : 0.00000000E+00 fermie : -3.68080047E-01 cartesian_stress_tensor: # hartree/bohr^3 - [ -4.02597615E-07, 0.00000000E+00, 0.00000000E+00, ] - [ 0.00000000E+00, 4.75261789E-06, 0.00000000E+00, ] - [ 0.00000000E+00, 0.00000000E+00, 4.75261789E-06, ] pressure_GPa: -8.9270E-02 xred : - [ -5.9938E-02, 0.0000E+00, 0.0000E+00, H] - [ 5.9938E-02, 0.0000E+00, 0.0000E+00, H] cartesian_forces: # hartree/bohr - [ -4.96425578E-03, -0.00000000E+00, -0.00000000E+00, ] - [ 4.96425578E-03, -0.00000000E+00, -0.00000000E+00, ] force_length_stats: {min: 4.96425578E-03, max: 4.96425578E-03, mean: 4.96425578E-03, } ... Integrated electronic density in atomic spheres: ------------------------------------------------ Atom Sphere_radius Integrated_density 1 2.00000 1.46016254 2 2.00000 1.46016254 ---OUTPUT----------------------------------------------------------------------- Cartesian coordinates (xcart) [bohr] -7.19259203377435E-01 0.00000000000000E+00 0.00000000000000E+00 7.19259203377435E-01 0.00000000000000E+00 0.00000000000000E+00 Reduced coordinates (xred) -5.99382669481196E-02 0.00000000000000E+00 0.00000000000000E+00 5.99382669481196E-02 0.00000000000000E+00 0.00000000000000E+00 Cartesian forces (fcart) [Ha/bohr]; max,rms= 4.96426E-03 2.86611E-03 (free atoms) -4.96425577745528E-03 -0.00000000000000E+00 -0.00000000000000E+00 4.96425577745528E-03 -0.00000000000000E+00 -0.00000000000000E+00 Reduced forces (fred) 5.95710693294634E-02 0.00000000000000E+00 0.00000000000000E+00 -5.95710693294634E-02 -0.00000000000000E+00 -0.00000000000000E+00 Total energy (etotal) [Ha]= -1.13290104758927E+00 Difference of energy with previous step (new-old): Absolute (Ha)=-4.60324E-04 Relative =-4.06406E-04 --- Iteration: ( 3/10) Internal Cycle: (1/1) -------------------------------------------------------------------------------- ---SELF-CONSISTENT-FIELD CONVERGENCE-------------------------------------------- --- !BeginCycle iteration_state: {dtset: 11, itime: 3, icycle: 1, } solver: {iscf: 7, nstep: 10, nline: 4, wfoptalg: 0, } tolerances: {toldff: 5.00E-05, } ... iter Etot(hartree) deltaE(h) residm vres2 diffor maxfor ETOT 1 -1.1329364993527 -1.133E+00 1.710E-11 1.107E-02 4.835E-03 1.295E-04 ETOT 2 -1.1329368709556 -3.716E-07 2.063E-13 4.015E-04 2.315E-04 3.610E-04 ETOT 3 -1.1329369175438 -4.659E-08 5.714E-11 1.412E-05 8.431E-05 4.453E-04 ETOT 4 -1.1329369190174 -1.474E-09 1.105E-11 1.469E-08 1.023E-05 4.351E-04 ETOT 5 -1.1329369190175 -1.337E-13 3.478E-16 2.187E-09 2.190E-07 4.348E-04 At SCF step 5, forces are converged : for the second time, max diff in force= 2.190E-07 < toldff= 5.000E-05 Cartesian components of stress tensor (hartree/bohr^3) sigma(1 1)= 3.38319928E-06 sigma(3 2)= 0.00000000E+00 sigma(2 2)= 4.71425618E-06 sigma(3 1)= 0.00000000E+00 sigma(3 3)= 4.71425618E-06 sigma(2 1)= 0.00000000E+00 --- !ResultsGS iteration_state: {dtset: 11, itime: 3, icycle: 1, } comment : Summary of ground state results lattice_vectors: - [ 12.0000000, 0.0000000, 0.0000000, ] - [ 0.0000000, 12.0000000, 0.0000000, ] - [ 0.0000000, 0.0000000, 12.0000000, ] lattice_lengths: [ 12.00000, 12.00000, 12.00000, ] lattice_angles: [ 90.000, 90.000, 90.000, ] # degrees, (23, 13, 12) lattice_volume: 1.7280000E+03 convergence: {deltae: -1.337E-13, res2: 2.187E-09, residm: 3.478E-16, diffor: 2.190E-07, } etotal : -1.13293692E+00 entropy : 0.00000000E+00 fermie : -3.66880383E-01 cartesian_stress_tensor: # hartree/bohr^3 - [ 3.38319928E-06, 0.00000000E+00, 0.00000000E+00, ] - [ 0.00000000E+00, 4.71425618E-06, 0.00000000E+00, ] - [ 0.00000000E+00, 0.00000000E+00, 4.71425618E-06, ] pressure_GPa: -1.2564E-01 xred : - [ -6.0496E-02, 0.0000E+00, 0.0000E+00, H] - [ 6.0496E-02, 0.0000E+00, 0.0000E+00, H] cartesian_forces: # hartree/bohr - [ -4.34840581E-04, -0.00000000E+00, -0.00000000E+00, ] - [ 4.34840581E-04, -0.00000000E+00, -0.00000000E+00, ] force_length_stats: {min: 4.34840581E-04, max: 4.34840581E-04, mean: 4.34840581E-04, } ... Integrated electronic density in atomic spheres: ------------------------------------------------ Atom Sphere_radius Integrated_density 1 2.00000 1.45324627 2 2.00000 1.45324627 ---OUTPUT----------------------------------------------------------------------- Cartesian coordinates (xcart) [bohr] -7.25947413387786E-01 0.00000000000000E+00 0.00000000000000E+00 7.25947413387786E-01 0.00000000000000E+00 0.00000000000000E+00 Reduced coordinates (xred) -6.04956177823155E-02 0.00000000000000E+00 0.00000000000000E+00 6.04956177823155E-02 0.00000000000000E+00 0.00000000000000E+00 Cartesian forces (fcart) [Ha/bohr]; max,rms= 4.34841E-04 2.51055E-04 (free atoms) -4.34840580812663E-04 -0.00000000000000E+00 -0.00000000000000E+00 4.34840580812663E-04 -0.00000000000000E+00 -0.00000000000000E+00 Reduced forces (fred) 5.21808696975196E-03 0.00000000000000E+00 0.00000000000000E+00 -5.21808696975196E-03 -0.00000000000000E+00 -0.00000000000000E+00 Total energy (etotal) [Ha]= -1.13293691901749E+00 Difference of energy with previous step (new-old): Absolute (Ha)=-3.58714E-05 Relative =-3.16628E-05 At Broyd/MD step 3, gradients are converged : max grad (force/stress) = 4.3484E-04 < tolmxf= 5.0000E-04 ha/bohr (free atoms) ================================================================================ ----iterations are completed or convergence reached---- Mean square residual over all n,k,spin= 34.778E-17; max= 34.778E-17 reduced coordinates (array xred) for 2 atoms -0.060495617782 0.000000000000 0.000000000000 0.060495617782 0.000000000000 0.000000000000 rms dE/dt= 3.0127E-03; max dE/dt= 5.2181E-03; dE/dt below (all hartree) 1 0.005218086970 0.000000000000 0.000000000000 2 -0.005218086970 0.000000000000 0.000000000000 cartesian coordinates (angstrom) at end: 1 -0.38415482579968 0.00000000000000 0.00000000000000 2 0.38415482579968 0.00000000000000 0.00000000000000 cartesian forces (hartree/bohr) at end: 1 -0.00043484058081 -0.00000000000000 -0.00000000000000 2 0.00043484058081 -0.00000000000000 -0.00000000000000 frms,max,avg= 2.5105533E-04 4.3484058E-04 0.000E+00 0.000E+00 0.000E+00 h/b cartesian forces (eV/Angstrom) at end: 1 -0.02236039982509 -0.00000000000000 -0.00000000000000 2 0.02236039982509 -0.00000000000000 -0.00000000000000 frms,max,avg= 1.2909783E-02 2.2360400E-02 0.000E+00 0.000E+00 0.000E+00 e/A length scales= 12.000000000000 12.000000000000 12.000000000000 bohr = 6.350126503080 6.350126503080 6.350126503080 angstroms prteigrs : about to open file tbase2_4o_DS11_EIG Fermi (or HOMO) energy (hartree) = -0.36688 Average Vxc (hartree)= -0.04722 Eigenvalues (hartree) for nkpt= 1 k points: kpt# 1, nband= 1, wtk= 1.00000, kpt= 0.0000 0.0000 0.0000 (reduced coord) -0.36688 --- !EnergyTerms iteration_state : {dtset: 11, itime: 3, icycle: 1, } comment : Components of total free energy in Hartree kinetic : 1.06135735371521E+00 hartree : 8.14007442739970E-01 xc : -6.41716190451920E-01 Ewald energy : 2.18482845469349E-01 psp_core : -1.11193990319380E-05 local_psp : -2.58505725109107E+00 non_local_psp : 0.00000000000000E+00 total_energy : -1.13293691901749E+00 total_energy_eV : -3.08287813925508E+01 band_energy : -7.33760765552470E-01 ... rms coord change= 1.2484E-03 atom, delta coord (reduced): 1 -0.002162284449 0.000000000000 0.000000000000 2 0.002162284449 0.000000000000 0.000000000000 Cartesian components of stress tensor (hartree/bohr^3) sigma(1 1)= 3.38319928E-06 sigma(3 2)= 0.00000000E+00 sigma(2 2)= 4.71425618E-06 sigma(3 1)= 0.00000000E+00 sigma(3 3)= 4.71425618E-06 sigma(2 1)= 0.00000000E+00 -Cartesian components of stress tensor (GPa) [Pressure= -1.2564E-01 GPa] - sigma(1 1)= 9.95371425E-02 sigma(3 2)= 0.00000000E+00 - sigma(2 2)= 1.38698182E-01 sigma(3 1)= 0.00000000E+00 - sigma(3 3)= 1.38698182E-01 sigma(2 1)= 0.00000000E+00 ================================================================================ == DATASET 12 ================================================================== - mpi_nproc: 1, omp_nthreads: -1 (-1 if OMP is not activated) --- !DatasetInfo iteration_state: {dtset: 12, } dimensions: {natom: 1, nkpt: 1, mband: 1, nsppol: 2, nspinor: 1, nspden: 2, mpw: 6759, } cutoff_energies: {ecut: 30.0, pawecutdg: -1.0, } electrons: {nelect: 1.00000000E+00, charge: 0.00000000E+00, occopt: 2.00000000E+00, tsmear: 1.00000000E-02, } meta: {optdriver: 0, ionmov: 0, optcell: 0, iscf: 7, paral_kgb: 0, } ... Exchange-correlation functional for the present dataset will be: LDA: new Teter (4/93) with spin-polarized option - ixc=1 Citation for XC functional: S. Goedecker, M. Teter, J. Huetter, PRB 54, 1703 (1996) Real(R)+Recip(G) space primitive vectors, cartesian coordinates (Bohr,Bohr^-1): R(1)= 12.0000000 0.0000000 0.0000000 G(1)= 0.0833333 0.0000000 0.0000000 R(2)= 0.0000000 12.0000000 0.0000000 G(2)= 0.0000000 0.0833333 0.0000000 R(3)= 0.0000000 0.0000000 12.0000000 G(3)= 0.0000000 0.0000000 0.0833333 Unit cell volume ucvol= 1.7280000E+03 bohr^3 Angles (23,13,12)= 9.00000000E+01 9.00000000E+01 9.00000000E+01 degrees getcut: wavevector= 0.0000 0.0000 0.0000 ngfft= 60 60 60 ecut(hartree)= 30.000 => boxcut(ratio)= 2.02789 -4.80358038E-03 ecore*ucvol(ha*bohr**3) -------------------------------------------------------------------------------- _setup2: Arith. and geom. avg. npw (full set) are 13517.000 13517.000 ================================================================================ --- !BeginCycle iteration_state: {dtset: 12, } solver: {iscf: 7, nstep: 10, nline: 4, wfoptalg: 0, } tolerances: {toldfe: 1.00E-06, } ... iter Etot(hartree) deltaE(h) residm vres2 ETOT 1 -0.47735047307878 -4.774E-01 6.905E-05 6.866E+01 ETOT 2 -0.47765217972552 -3.017E-04 3.555E-11 2.284E-01 ETOT 3 -0.47765320174752 -1.022E-06 6.783E-09 3.454E-04 ETOT 4 -0.47765320720361 -5.456E-09 8.347E-11 1.477E-06 ETOT 5 -0.47765320721318 -9.577E-12 9.542E-14 8.495E-09 At SCF step 5, etot is converged : for the second time, diff in etot= 9.577E-12 < toldfe= 1.000E-06 Cartesian components of stress tensor (hartree/bohr^3) sigma(1 1)= 2.39435055E-06 sigma(3 2)= 0.00000000E+00 sigma(2 2)= 2.39435055E-06 sigma(3 1)= 0.00000000E+00 sigma(3 3)= 2.39435055E-06 sigma(2 1)= 0.00000000E+00 --- !ResultsGS iteration_state: {dtset: 12, } comment : Summary of ground state results lattice_vectors: - [ 12.0000000, 0.0000000, 0.0000000, ] - [ 0.0000000, 12.0000000, 0.0000000, ] - [ 0.0000000, 0.0000000, 12.0000000, ] lattice_lengths: [ 12.00000, 12.00000, 12.00000, ] lattice_angles: [ 90.000, 90.000, 90.000, ] # degrees, (23, 13, 12) lattice_volume: 1.7280000E+03 convergence: {deltae: -9.577E-12, res2: 8.495E-09, residm: 9.542E-14, diffor: null, } etotal : -4.77653207E-01 entropy : 0.00000000E+00 fermie : -2.66036564E-01 cartesian_stress_tensor: # hartree/bohr^3 - [ 2.39435055E-06, 0.00000000E+00, 0.00000000E+00, ] - [ 0.00000000E+00, 2.39435055E-06, 0.00000000E+00, ] - [ 0.00000000E+00, 0.00000000E+00, 2.39435055E-06, ] pressure_GPa: -7.0444E-02 xred : - [ 0.0000E+00, 0.0000E+00, 0.0000E+00, H] cartesian_forces: # hartree/bohr - [ -0.00000000E+00, -0.00000000E+00, -0.00000000E+00, ] force_length_stats: {min: 0.00000000E+00, max: 0.00000000E+00, mean: 0.00000000E+00, } ... Integrated electronic and magnetization densities in atomic spheres: --------------------------------------------------------------------- Radius=ratsph(iatom), smearing ratsm= 0.0000. Diff(up-dn)=approximate z local magnetic moment. Atom Radius up_density dn_density Total(up+dn) Diff(up-dn) 1 2.00000 0.722835 0.000000 0.722835 0.722835 --------------------------------------------------------------------- Sum: 0.722835 0.000000 0.722835 0.722835 Total magnetization (from the atomic spheres): 0.722835 Total magnetization (exact up - dn): 1.000000 ================================================================================ ----iterations are completed or convergence reached---- Mean square residual over all n,k,spin= 67.086E-15; max= 95.421E-15 reduced coordinates (array xred) for 1 atoms 0.000000000000 0.000000000000 0.000000000000 rms dE/dt= 0.0000E+00; max dE/dt= 0.0000E+00; dE/dt below (all hartree) 1 0.000000000000 0.000000000000 0.000000000000 cartesian coordinates (angstrom) at end: 1 0.00000000000000 0.00000000000000 0.00000000000000 cartesian forces (hartree/bohr) at end: 1 -0.00000000000000 -0.00000000000000 -0.00000000000000 frms,max,avg= 0.0000000E+00 0.0000000E+00 0.000E+00 0.000E+00 0.000E+00 h/b cartesian forces (eV/Angstrom) at end: 1 -0.00000000000000 -0.00000000000000 -0.00000000000000 frms,max,avg= 0.0000000E+00 0.0000000E+00 0.000E+00 0.000E+00 0.000E+00 e/A length scales= 12.000000000000 12.000000000000 12.000000000000 bohr = 6.350126503080 6.350126503080 6.350126503080 angstroms prteigrs : about to open file tbase2_4o_DS12_EIG Fermi (or HOMO) energy (hartree) = -0.26604 Average Vxc (hartree)= -0.04490 Eigenvalues (hartree) for nkpt= 1 k points, SPIN UP: kpt# 1, nband= 1, wtk= 1.00000, kpt= 0.0000 0.0000 0.0000 (reduced coord) -0.26604 Eigenvalues (hartree) for nkpt= 1 k points, SPIN DOWN: kpt# 1, nband= 1, wtk= 1.00000, kpt= 0.0000 0.0000 0.0000 (reduced coord) -0.10016 --- !EnergyTerms iteration_state : {dtset: 12, } comment : Components of total free energy in Hartree kinetic : 4.54200528773610E-01 hartree : 1.81130249718783E-01 xc : -2.75288191622281E-01 Ewald energy : -1.18220728311694E-01 psp_core : -2.77984975798451E-06 local_psp : -7.19472285921843E-01 non_local_psp : 0.00000000000000E+00 total_energy : -4.77653207213182E-01 total_energy_eV : -1.29976047734380E+01 band_energy : -2.66036564299337E-01 ... Cartesian components of stress tensor (hartree/bohr^3) sigma(1 1)= 2.39435055E-06 sigma(3 2)= 0.00000000E+00 sigma(2 2)= 2.39435055E-06 sigma(3 1)= 0.00000000E+00 sigma(3 3)= 2.39435055E-06 sigma(2 1)= 0.00000000E+00 -Cartesian components of stress tensor (GPa) [Pressure= -7.0444E-02 GPa] - sigma(1 1)= 7.04442134E-02 sigma(3 2)= 0.00000000E+00 - sigma(2 2)= 7.04442134E-02 sigma(3 1)= 0.00000000E+00 - sigma(3 3)= 7.04442134E-02 sigma(2 1)= 0.00000000E+00 == END DATASET(S) ============================================================== ================================================================================ -outvars: echo values of variables after computation -------- acell 1.2000000000E+01 1.2000000000E+01 1.2000000000E+01 Bohr amu 1.00794000E+00 bs_loband11 0 bs_loband12 0 0 diemac 1.00000000E+00 ecut 3.00000000E+01 Hartree etotal11 -1.1329369190E+00 etotal12 -4.7765320721E-01 fcart11 -4.3484058081E-04 -0.0000000000E+00 -0.0000000000E+00 4.3484058081E-04 -0.0000000000E+00 -0.0000000000E+00 fcart12 -0.0000000000E+00 -0.0000000000E+00 -0.0000000000E+00 - fftalg 312 ionmov11 2 ionmov12 0 istwfk 2 jdtset 11 12 kptopt 0 P mkmem 1 natom11 2 natom12 1 nband11 1 nband12 1 1 ndtset 2 ngfft 60 60 60 nkpt 1 nspden11 1 nspden12 2 nsppol11 1 nsppol12 2 nstep 10 nsym11 16 nsym12 48 ntime11 10 ntime12 1 ntypat 1 occ11 2.000000 occ12 1.000000 0.000000 occopt11 1 occopt12 2 optforces11 1 optforces12 2 spgroup11 123 spgroup12 221 spinat12 0.0000000000E+00 0.0000000000E+00 1.0000000000E+00 strten11 3.3831992765E-06 4.7142561810E-06 4.7142561810E-06 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 strten12 2.3943505506E-06 2.3943505506E-06 2.3943505506E-06 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 symafm11 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 symafm12 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 symrel11 1 0 0 0 1 0 0 0 1 -1 0 0 0 -1 0 0 0 -1 -1 0 0 0 1 0 0 0 -1 1 0 0 0 -1 0 0 0 1 -1 0 0 0 -1 0 0 0 1 1 0 0 0 1 0 0 0 -1 1 0 0 0 -1 0 0 0 -1 -1 0 0 0 1 0 0 0 1 1 0 0 0 0 1 0 1 0 -1 0 0 0 0 -1 0 -1 0 -1 0 0 0 0 1 0 -1 0 1 0 0 0 0 -1 0 1 0 -1 0 0 0 0 -1 0 1 0 1 0 0 0 0 1 0 -1 0 1 0 0 0 0 -1 0 -1 0 -1 0 0 0 0 1 0 1 0 symrel12 1 0 0 0 1 0 0 0 1 -1 0 0 0 -1 0 0 0 -1 -1 0 0 0 1 0 0 0 -1 1 0 0 0 -1 0 0 0 1 -1 0 0 0 -1 0 0 0 1 1 0 0 0 1 0 0 0 -1 1 0 0 0 -1 0 0 0 -1 -1 0 0 0 1 0 0 0 1 0 1 0 1 0 0 0 0 1 0 -1 0 -1 0 0 0 0 -1 0 -1 0 1 0 0 0 0 -1 0 1 0 -1 0 0 0 0 1 0 -1 0 -1 0 0 0 0 1 0 1 0 1 0 0 0 0 -1 0 1 0 -1 0 0 0 0 -1 0 -1 0 1 0 0 0 0 1 0 0 1 1 0 0 0 1 0 0 0 -1 -1 0 0 0 -1 0 0 0 -1 1 0 0 0 -1 0 0 0 1 -1 0 0 0 1 0 0 0 -1 -1 0 0 0 1 0 0 0 1 1 0 0 0 -1 0 0 0 1 -1 0 0 0 -1 0 0 0 -1 1 0 0 0 1 0 1 0 0 0 0 1 0 1 0 -1 0 0 0 0 -1 0 -1 0 -1 0 0 0 0 1 0 -1 0 1 0 0 0 0 -1 0 1 0 -1 0 0 0 0 -1 0 1 0 1 0 0 0 0 1 0 -1 0 1 0 0 0 0 -1 0 -1 0 -1 0 0 0 0 1 0 1 0 0 1 0 0 0 1 1 0 0 0 -1 0 0 0 -1 -1 0 0 0 -1 0 0 0 1 -1 0 0 0 1 0 0 0 -1 1 0 0 0 -1 0 0 0 -1 1 0 0 0 1 0 0 0 1 -1 0 0 0 1 0 0 0 -1 -1 0 0 0 -1 0 0 0 1 1 0 0 0 0 1 0 1 0 1 0 0 0 0 -1 0 -1 0 -1 0 0 0 0 -1 0 1 0 -1 0 0 0 0 1 0 -1 0 1 0 0 0 0 -1 0 -1 0 1 0 0 0 0 1 0 1 0 -1 0 0 0 0 1 0 -1 0 -1 0 0 0 0 -1 0 1 0 1 0 0 tnons11 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 tnons12 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 0.0000000 toldfe11 0.00000000E+00 Hartree toldfe12 1.00000000E-06 Hartree toldff11 5.00000000E-05 toldff12 0.00000000E+00 tolmxf11 5.00000000E-04 tolmxf12 5.00000000E-05 typat11 1 1 typat12 1 xangst11 -3.8415482580E-01 0.0000000000E+00 0.0000000000E+00 3.8415482580E-01 0.0000000000E+00 0.0000000000E+00 xangst12 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 xcart11 -7.2594741339E-01 0.0000000000E+00 0.0000000000E+00 7.2594741339E-01 0.0000000000E+00 0.0000000000E+00 xcart12 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 xred11 -6.0495617782E-02 0.0000000000E+00 0.0000000000E+00 6.0495617782E-02 0.0000000000E+00 0.0000000000E+00 xred12 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00 znucl 1.00000 ================================================================================ - Timing analysis has been suppressed with timopt=0 ================================================================================ Suggested references for the acknowledgment of ABINIT usage. The users of ABINIT have little formal obligations with respect to the ABINIT group (those specified in the GNU General Public License, http://www.gnu.org/copyleft/gpl.txt). However, it is common practice in the scientific literature, to acknowledge the efforts of people that have made the research possible. In this spirit, please find below suggested citations of work written by ABINIT developers, corresponding to implementations inside of ABINIT that you have used in the present run. Note also that it will be of great value to readers of publications presenting these results, to read papers enabling them to understand the theoretical formalism and details of the ABINIT implementation. For information on why they are suggested, see also https://docs.abinit.org/theory/acknowledgments. - - [1] The Abinit project: Impact, environment and recent developments. - Computer Phys. Comm. 248, 107042 (2020). - X.Gonze, B. Amadon, G. Antonius, F.Arnardi, L.Baguet, J.-M.Beuken, - J.Bieder, F.Bottin, J.Bouchet, E.Bousquet, N.Brouwer, F.Bruneval, - G.Brunin, T.Cavignac, J.-B. Charraud, Wei Chen, M.Cote, S.Cottenier, - J.Denier, G.Geneste, Ph.Ghosez, M.Giantomassi, Y.Gillet, O.Gingras, - D.R.Hamann, G.Hautier, Xu He, N.Helbig, N.Holzwarth, Y.Jia, F.Jollet, - W.Lafargue-Dit-Hauret, K.Lejaeghere, M.A.L.Marques, A.Martin, C.Martins, - H.P.C. Miranda, F.Naccarato, K. Persson, G.Petretto, V.Planes, Y.Pouillon, - S.Prokhorenko, F.Ricci, G.-M.Rignanese, A.H.Romero, M.M.Schmitt, M.Torrent, - M.J.van Setten, B.Van Troeye, M.J.Verstraete, G.Zerah and J.W.Zwanzig - Comment: the fifth generic paper describing the ABINIT project. - Note that a version of this paper, that is not formatted for Computer Phys. Comm. - is available at https://www.abinit.org/sites/default/files/ABINIT20.pdf . - The licence allows the authors to put it on the Web. - DOI and bibtex: see https://docs.abinit.org/theory/bibliography/#gonze2020 - - [2] Recent developments in the ABINIT software package. - Computer Phys. Comm. 205, 106 (2016). - X.Gonze, F.Jollet, F.Abreu Araujo, D.Adams, B.Amadon, T.Applencourt, - C.Audouze, J.-M.Beuken, J.Bieder, A.Bokhanchuk, E.Bousquet, F.Bruneval - D.Caliste, M.Cote, F.Dahm, F.Da Pieve, M.Delaveau, M.Di Gennaro, - B.Dorado, C.Espejo, G.Geneste, L.Genovese, A.Gerossier, M.Giantomassi, - Y.Gillet, D.R.Hamann, L.He, G.Jomard, J.Laflamme Janssen, S.Le Roux, - A.Levitt, A.Lherbier, F.Liu, I.Lukacevic, A.Martin, C.Martins, - M.J.T.Oliveira, S.Ponce, Y.Pouillon, T.Rangel, G.-M.Rignanese, - A.H.Romero, B.Rousseau, O.Rubel, A.A.Shukri, M.Stankovski, M.Torrent, - M.J.Van Setten, B.Van Troeye, M.J.Verstraete, D.Waroquier, J.Wiktor, - B.Xu, A.Zhou, J.W.Zwanziger. - Comment: the fourth generic paper describing the ABINIT project. - Note that a version of this paper, that is not formatted for Computer Phys. Comm. - is available at https://www.abinit.org/sites/default/files/ABINIT16.pdf . - The licence allows the authors to put it on the Web. - DOI and bibtex: see https://docs.abinit.org/theory/bibliography/#gonze2016 - - [3] ABINIT: First-principles approach of materials and nanosystem properties. - Computer Phys. Comm. 180, 2582-2615 (2009). - X. Gonze, B. Amadon, P.-M. Anglade, J.-M. Beuken, F. Bottin, P. Boulanger, F. Bruneval, - D. Caliste, R. Caracas, M. Cote, T. Deutsch, L. Genovese, Ph. Ghosez, M. Giantomassi - S. Goedecker, D.R. Hamann, P. Hermet, F. Jollet, G. Jomard, S. Leroux, M. Mancini, S. Mazevet, - M.J.T. Oliveira, G. Onida, Y. Pouillon, T. Rangel, G.-M. Rignanese, D. Sangalli, R. Shaltaf, - M. Torrent, M.J. Verstraete, G. Zerah, J.W. Zwanziger - Comment: the third generic paper describing the ABINIT project. - Note that a version of this paper, that is not formatted for Computer Phys. Comm. - is available at https://www.abinit.org/sites/default/files/ABINIT_CPC_v10.pdf . - The licence allows the authors to put it on the Web. - DOI and bibtex: see https://docs.abinit.org/theory/bibliography/#gonze2009 - - And optionally: - - [4] A brief introduction to the ABINIT software package. - Z. Kristallogr. 220, 558-562 (2005). - X. Gonze, G.-M. Rignanese, M. Verstraete, J.-M. Beuken, Y. Pouillon, R. Caracas, F. Jollet, - M. Torrent, G. Zerah, M. Mikami, Ph. Ghosez, M. Veithen, J.-Y. Raty, V. Olevano, F. Bruneval, - L. Reining, R. Godby, G. Onida, D.R. Hamann, and D.C. Allan. - Comment: the second generic paper describing the ABINIT project. Note that this paper - should be cited especially if you are using the GW part of ABINIT, as several authors - of this part are not in the list of authors of the first or third paper. - The .pdf of the latter paper is available at https://www.abinit.org/sites/default/files/zfk_0505-06_558-562.pdf. - Note that it should not redistributed (Copyright by Oldenburg Wissenschaftverlag, - the licence allows the authors to put it on the Web). - DOI and bibtex: see https://docs.abinit.org/theory/bibliography/#gonze2005 - - Proc. 0 individual time (sec): cpu= 3.6 wall= 3.8 ================================================================================ Calculation completed. .Delivered 5 WARNINGs and 8 COMMENTs to log file. +Overall time at end (sec) : cpu= 3.6 wall= 3.8
Since we are doing the calculation at a single (ecut, acell) pair, the total CPU time is not as much as for the previous determinations of optimal values through series calculations. However, the memory needs have still increased a bit.
The output data are:
etotal11 -1.1329369190E+00 etotal12 -4.7765320721E-01 xcart11 -7.2594741339E-01 0.0000000000E+00 0.0000000000E+00 7.2594741339E-01 0.0000000000E+00 0.0000000000E+00 xcart12 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00
- The corresponding atomisation energy is
0.1776 Ha = 4.833 eV
- The interatomic distance is 1.452 Bohr.
- These are our final data for the local (spin) density approximation.
We have used ixc = 1. Other expressions for the local (spin) density approximation [2, 3 … 7] are possible. The values 1, 2, 3 and 7 should give about the same results, since they all start from the XC energy of the homogeneous electron gas, as determined by Quantum Monte Carlo calculations. Other possibilities (ixc = 4, 5, 6) are older local density functionals, that could not rely on these data.
6 The use of the Generalized Gradient Approximation¶
We will use the Perdew-Burke-Ernzerhof functional proposed in [Perdew1996]
In principle, for GGA, one should use another pseudopotential than for LDA. However, for the special case of Hydrogen, and in general pseudopotentials with a very small core (including only the 1s orbital), pseudopotentials issued from the LDA and from the GGA are very similar. So, we will not change our pseudopotential. This will save us lot of time, as we should not redo an ecut convergence test
Important
ecut is often characteristic of the pseudopotentials that are used in a calculation.
Independently of the pseudopotential, an acell convergence test should not be done again, since the vacuum is treated similarly in LDA or GGA.
So, our final values within GGA will be easily obtained by setting ixc to 11 in tbase2_4.in.
# H2 molecule in a big box # # Like t24.in, except that GGA is used instead of LDA. # ndtset 2 udtset 1 2 # Use the Perdew-Burke-Ernzerhof GGA ixc 11 #Definition of the unit cell and ecut, #for which one will have to make a convergence study acell:? 12 12 12 acell+? 2 2 2 ecut 30 #First dataset : find the optimal bond length of H2, and associated total energy natom?1 2 # There are two atoms ionmov?1 2 # Use the modified Broyden algorithm ntime?1 10 # Maximum number of Broyden "timesteps" tolmxf?1 5.0d-4 # Stopping criterion for the geometry optimization : when # the residual forces are less than tolmxf, the Broyden # algorithm can stop xcart?1 -0.7 0.0 0.0 # The starting values of the 0.7 0.0 0.0 # atomic coordinates toldff?1 5.0d-5 # Will stop the SCF cycle when, twice in a row, # the difference between two consecutive evaluations of # forces differ by less than toldff (in Hartree/Bohr) nband?1 1 # Just one band #Second dataset : get the total energy of the isolated atom natom?2 1 # There is one atom nsppol?2 2 # Spin-polarized calculation occopt?2 2 # Allow occupation numbers to be set by hand nband?2 1 1 # Number of bands for spin up and spin down occ?2 1.0 0.0 # Occupation numbers for spin up state and spin down state. toldfe?2 1.0d-6 # Will stop the SCF cycles when, twice in a row, # the difference between two consecutive evaluations # of total energy differ by less than toldfe (in Hartree) # This value is way too large for most realistic studies of materials xcart?2 0.0 0.0 0.0 # The atom is located at the origin spinat?2 0.0 0.0 1.0 # Initialisation of spin #rprim 1 0 0 0 1 0 0 0 1 # This line, defining orthogonal primitive vectors, # is commented, because it is precisely the default value of rprim #Definition of the atom types ntypat 1 # There is only one type of atom znucl 1 # The keyword "znucl" refers to the atomic number of the # possible type(s) of atom. The pseudopotential(s) # mentioned in the "files" file must correspond # to the type(s) of atom. Here, the only type is Hydrogen. #Definition of the atoms typat 1 1 # For the first dataset, both numbers will be read, # while for the second dataset, only one number will be read #Definition of the k-point grid kptopt 0 # Enter the k points manually nkpt 1 # Only one k point is needed for isolated system, # taken by default to be 0.0 0.0 0.0 #Definition of the SCF procedure nstep 10 # Maximal number of SCF cycles #toldfe is no more defined, as toldff is used above... diemac 2.0 # Although this is not mandatory, it is worth to # precondition the SCF cycle. The model dielectric # function used as the standard preconditioner # is described in the "dielng" input variable section. # Here, we follow the prescriptions for molecules # in a big box #%%<BEGIN TEST_INFO> #%% [setup] #%% executable = abinit #%% [files] #%% files_to_test = #%% tbase2_5.out, tolnlines= 0, tolabs= 3.681e-09, tolrel= 2.000e-04, fld_options = -medium #%% psp_files = 01h.pspgth #%% [paral_info] #%% max_nprocs = 1 #%% [extra_info] #%% authors = Unknown #%% keywords = #%% description = #%% H2 molecule in a big box #%% Like t24.in, except that GGA is used instead of LDA. #%%<END TEST_INFO>
etotal11 -1.1621428376E+00 etotal12 -4.9869631917E-01 xcart11 -7.1190611804E-01 0.0000000000E+00 0.0000000000E+00 7.1190611804E-01 0.0000000000E+00 0.0000000000E+00 xcart12 0.0000000000E+00 0.0000000000E+00 0.0000000000E+00
- The corresponding atomisation energy is 0.1648 Ha = 4.483 eV
- The interatomic distance is 1.424 Bohr.
- These are our final data for the generalized gradient approximation.
Once more, here are the experimental data:
- bond length: 1.401 Bohr
- atomisation energy: 4.747 eV
In GGA, we are within 2% of the experimental bond length, but 5% of the experimental atomisation energy. In LDA, we were within 4% of the experimental bond length, and within 2% of the experimental atomisation energy.
Important
Do not forget that the typical accuracy of LDA and GGA varies with the class of materials studied.